xref: /libCEED/rust/libceed-sys/c-src/backends/hip-gen/ceed-hip-gen-operator-build.cpp (revision f815fac990b20019e227e6950cc74a96439f9eba)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
37d8d0e25Snbeams //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
57d8d0e25Snbeams //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d576824SJeremy L Thompson 
87d8d0e25Snbeams #define CEED_DEBUG_COLOR 12
97d8d0e25Snbeams 
1049aac155SJeremy L Thompson #include <ceed.h>
11ec3da8bcSJed Brown #include <ceed/backend.h>
129e201c85SYohann #include <ceed/jit-tools.h>
132b730f8bSJeremy L Thompson 
147d8d0e25Snbeams #include <iostream>
157d8d0e25Snbeams #include <sstream>
162b730f8bSJeremy L Thompson #include <string>
172b730f8bSJeremy L Thompson 
180d0321e0SJeremy L Thompson #include "../hip-ref/ceed-hip-ref.h"
197d8d0e25Snbeams #include "../hip-shared/ceed-hip-shared.h"
20b2165e7aSSebastian Grimberg #include "../hip/ceed-hip-common.h"
217d8d0e25Snbeams #include "../hip/ceed-hip-compile.h"
222b730f8bSJeremy L Thompson #include "ceed-hip-gen.h"
23b3e1519bSnbeams 
24b3e1519bSnbeams //------------------------------------------------------------------------------
25b3e1519bSnbeams // Calculate the block size used for launching the operator kernel
26b3e1519bSnbeams //------------------------------------------------------------------------------
272b730f8bSJeremy L Thompson extern "C" int BlockGridCalculate_Hip_gen(const CeedInt dim, const CeedInt num_elem, const CeedInt P_1d, const CeedInt Q_1d, CeedInt *block_sizes) {
289e201c85SYohann   const CeedInt thread1d = CeedIntMax(Q_1d, P_1d);
29b3e1519bSnbeams   if (dim == 1) {
309e201c85SYohann     CeedInt elems_per_block = 64 * thread1d > 256 ? 256 / thread1d : 64;
31b7453713SJeremy L Thompson 
329e201c85SYohann     elems_per_block = elems_per_block > 0 ? elems_per_block : 1;
33b3e1519bSnbeams     block_sizes[0]  = thread1d;
34b3e1519bSnbeams     block_sizes[1]  = 1;
359e201c85SYohann     block_sizes[2]  = elems_per_block;
36b3e1519bSnbeams   } else if (dim == 2) {
379e201c85SYohann     const CeedInt elems_per_block = thread1d < 4 ? 16 : 2;
38b7453713SJeremy L Thompson 
39b3e1519bSnbeams     block_sizes[0] = thread1d;
40b3e1519bSnbeams     block_sizes[1] = thread1d;
419e201c85SYohann     block_sizes[2] = elems_per_block;
42b3e1519bSnbeams   } else if (dim == 3) {
439e201c85SYohann     const CeedInt elems_per_block = thread1d < 6 ? 4 : (thread1d < 8 ? 2 : 1);
44b7453713SJeremy L Thompson 
45b3e1519bSnbeams     block_sizes[0] = thread1d;
46b3e1519bSnbeams     block_sizes[1] = thread1d;
479e201c85SYohann     block_sizes[2] = elems_per_block;
48b3e1519bSnbeams   }
49b3e1519bSnbeams   return CEED_ERROR_SUCCESS;
50b3e1519bSnbeams }
51b3e1519bSnbeams 
527d8d0e25Snbeams //------------------------------------------------------------------------------
534b3e95d5SJeremy L Thompson // Determine type of operator
544b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
554b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Hip_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
564b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
574b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor,
584b3e95d5SJeremy L Thompson                                                bool *use_3d_slices) {
594b3e95d5SJeremy L Thompson   // Find dim, P_1d, Q_1d
604b3e95d5SJeremy L Thompson   *max_P_1d  = 0;
614b3e95d5SJeremy L Thompson   *Q_1d      = 0;
624b3e95d5SJeremy L Thompson   *dim       = 0;
634b3e95d5SJeremy L Thompson   *is_tensor = true;
644b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
654b3e95d5SJeremy L Thompson     CeedBasis basis;
664b3e95d5SJeremy L Thompson 
674b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
684b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
694b3e95d5SJeremy L Thompson       bool    is_field_tensor;
704b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
714b3e95d5SJeremy L Thompson 
724b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
734b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
744b3e95d5SJeremy L Thompson       CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
754b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
764b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
774b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
784b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
794b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
804b3e95d5SJeremy L Thompson       *dim = field_dim;
814b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
824b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
834b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
844b3e95d5SJeremy L Thompson     }
854b3e95d5SJeremy L Thompson   }
864b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
874b3e95d5SJeremy L Thompson     CeedBasis basis;
884b3e95d5SJeremy L Thompson 
894b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
904b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
914b3e95d5SJeremy L Thompson       bool    is_field_tensor;
924b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
934b3e95d5SJeremy L Thompson 
944b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
954b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
964b3e95d5SJeremy L Thompson       CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
974b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
984b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
994b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
1004b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
1014b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
1024b3e95d5SJeremy L Thompson       *dim = field_dim;
1034b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
1044b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
1054b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
1064b3e95d5SJeremy L Thompson     }
1074b3e95d5SJeremy L Thompson   }
1084b3e95d5SJeremy L Thompson 
1094b3e95d5SJeremy L Thompson   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1104b3e95d5SJeremy L Thompson   *use_3d_slices = false;
1114b3e95d5SJeremy L Thompson   if (*dim == 3) {
1124b3e95d5SJeremy L Thompson     bool was_grad_found = false;
1134b3e95d5SJeremy L Thompson 
1144b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1154b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1164b3e95d5SJeremy L Thompson 
1174b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1184b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1194b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
1204b3e95d5SJeremy L Thompson         CeedBasis             basis;
1214b3e95d5SJeremy L Thompson 
1224b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1234b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1244b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1254b3e95d5SJeremy L Thompson         was_grad_found = true;
1264b3e95d5SJeremy L Thompson       }
1274b3e95d5SJeremy L Thompson     }
1284b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1294b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1304b3e95d5SJeremy L Thompson 
1314b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1324b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1334b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
1344b3e95d5SJeremy L Thompson         CeedBasis             basis;
1354b3e95d5SJeremy L Thompson 
1364b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1374b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1384b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1394b3e95d5SJeremy L Thompson         was_grad_found = true;
1404b3e95d5SJeremy L Thompson       }
1414b3e95d5SJeremy L Thompson     }
1424b3e95d5SJeremy L Thompson   }
1434b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1444b3e95d5SJeremy L Thompson }
1454b3e95d5SJeremy L Thompson 
1464b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1474b3e95d5SJeremy L Thompson // Setup fields
1484b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1494b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedOperatorField op_field,
1504b3e95d5SJeremy L Thompson                                                     CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool use_3d_slices) {
1514b3e95d5SJeremy L Thompson   std::string           var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
1524b3e95d5SJeremy L Thompson   std::string           P_name = "P_1d" + var_suffix, Q_name = "Q_1d";
1534b3e95d5SJeremy L Thompson   std::string           option_name = (is_input ? "inputs" : "outputs");
1544b3e95d5SJeremy L Thompson   CeedEvalMode          eval_mode   = CEED_EVAL_NONE;
1554b3e95d5SJeremy L Thompson   CeedInt               elem_size = 0, num_comp = 0, P_1d = 0;
1564b3e95d5SJeremy L Thompson   CeedElemRestriction   elem_rstr;
1574b3e95d5SJeremy L Thompson   CeedBasis_Hip_shared *basis_data;
1584b3e95d5SJeremy L Thompson   CeedBasis             basis;
1594b3e95d5SJeremy L Thompson 
1604b3e95d5SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n";
1614b3e95d5SJeremy L Thompson 
1624b3e95d5SJeremy L Thompson   // Get field data
1634b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
1644b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
1654b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1664b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1674b3e95d5SJeremy L Thompson   }
1684b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
1694b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
1704b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1714b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1724b3e95d5SJeremy L Thompson   }
1734b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
1744b3e95d5SJeremy L Thompson 
1754b3e95d5SJeremy L Thompson   // Set field constants
1764b3e95d5SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1774b3e95d5SJeremy L Thompson     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
1784b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
1794b3e95d5SJeremy L Thompson   }
1804b3e95d5SJeremy L Thompson 
1814b3e95d5SJeremy L Thompson   // Load basis data
1824b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
1834b3e95d5SJeremy L Thompson   switch (eval_mode) {
1844b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
1854b3e95d5SJeremy L Thompson       break;
1864b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
1874b3e95d5SJeremy L Thompson       if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
1884b3e95d5SJeremy L Thompson       else data->B.outputs[i] = basis_data->d_interp_1d;
1894b3e95d5SJeremy L Thompson       code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n";
190*f815fac9SJeremy L Thompson       code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
1914b3e95d5SJeremy L Thompson       break;
1924b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
1934b3e95d5SJeremy L Thompson       if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
1944b3e95d5SJeremy L Thompson       else data->B.outputs[i] = basis_data->d_interp_1d;
1954b3e95d5SJeremy L Thompson       code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n";
196*f815fac9SJeremy L Thompson       code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
1974b3e95d5SJeremy L Thompson       if (use_3d_slices) {
1984b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
1994b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
2004b3e95d5SJeremy L Thompson         code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n";
201*f815fac9SJeremy L Thompson         code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2024b3e95d5SJeremy L Thompson       } else {
2034b3e95d5SJeremy L Thompson         bool has_collo_grad = basis_data->d_collo_grad_1d;
2044b3e95d5SJeremy L Thompson 
2054b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2064b3e95d5SJeremy L Thompson         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2074b3e95d5SJeremy L Thompson         if (has_collo_grad) {
2084b3e95d5SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n";
209*f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2104b3e95d5SJeremy L Thompson         } else {
2114b3e95d5SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * P_1d << "];\n";
212*f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2134b3e95d5SJeremy L Thompson         }
2144b3e95d5SJeremy L Thompson       }
2154b3e95d5SJeremy L Thompson       break;
2164b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
2174b3e95d5SJeremy L Thompson       break;  // No action
2184b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
2194b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
2204b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
2214b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
2224b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
2234b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
2244b3e95d5SJeremy L Thompson   }
2254b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2264b3e95d5SJeremy L Thompson }
2274b3e95d5SJeremy L Thompson 
2284b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2294b3e95d5SJeremy L Thompson // Restriction
2304b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2314b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
232e93651e5SJeremy L Thompson                                                       CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
233e93651e5SJeremy L Thompson                                                       CeedInt Q_1d, bool is_input, bool use_3d_slices) {
2344b3e95d5SJeremy L Thompson   std::string              var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
2354b3e95d5SJeremy L Thompson   std::string              P_name     = "P_1d" + var_suffix;
2364b3e95d5SJeremy L Thompson   CeedEvalMode             eval_mode  = CEED_EVAL_NONE;
2374b3e95d5SJeremy L Thompson   CeedInt                  elem_size = 0, num_comp = 0, P_1d = 0;
2384b3e95d5SJeremy L Thompson   CeedSize                 l_size;
239*f815fac9SJeremy L Thompson   CeedRestrictionType      rstr_type = CEED_RESTRICTION_STANDARD;
2404b3e95d5SJeremy L Thompson   CeedElemRestriction_Hip *rstr_data;
2414b3e95d5SJeremy L Thompson   CeedElemRestriction      elem_rstr;
2424b3e95d5SJeremy L Thompson   CeedBasis                basis;
2434b3e95d5SJeremy L Thompson 
2444b3e95d5SJeremy L Thompson   // Get field data
2454b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
2464b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
247*f815fac9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
2484b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
2494b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
2504b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
2514b3e95d5SJeremy L Thompson   }
2524b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
2534b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
2544b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2554b3e95d5SJeremy L Thompson   }
2564b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
2574b3e95d5SJeremy L Thompson 
2584b3e95d5SJeremy L Thompson   // Restriction
2594b3e95d5SJeremy L Thompson   if (is_input) {
2604b3e95d5SJeremy L Thompson     // Input
261e93651e5SJeremy L Thompson     if (field_input_buffer[i] != i) {
262e93651e5SJeremy L Thompson       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
263e93651e5SJeremy L Thompson 
264e93651e5SJeremy L Thompson       // Restriction was already done for previous input
265e93651e5SJeremy L Thompson       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
266e93651e5SJeremy L Thompson     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices)) {
267e93651e5SJeremy L Thompson       if (eval_mode == CEED_EVAL_NONE) {
268e93651e5SJeremy L Thompson         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
2694b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
270e93651e5SJeremy L Thompson       } else {
271e93651e5SJeremy L Thompson         // Otherwise we're using the scratch space
272e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
273e93651e5SJeremy L Thompson       }
274*f815fac9SJeremy L Thompson       switch (rstr_type) {
275*f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STANDARD: {
2764b3e95d5SJeremy L Thompson           CeedInt comp_stride;
2774b3e95d5SJeremy L Thompson 
2784b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
2794b3e95d5SJeremy L Thompson           code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
2804b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
2814b3e95d5SJeremy L Thompson           code << "    // CompStride: " << comp_stride << "\n";
2824b3e95d5SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
283*f815fac9SJeremy L Thompson           code << "    ReadLVecStandard" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size"
284*f815fac9SJeremy L Thompson                << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
285*f815fac9SJeremy L Thompson           break;
286*f815fac9SJeremy L Thompson         }
287*f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STRIDED: {
2884b3e95d5SJeremy L Thompson           bool    has_backend_strides;
2894b3e95d5SJeremy L Thompson           CeedInt num_elem;
2904b3e95d5SJeremy L Thompson 
2914b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
2924b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
2934b3e95d5SJeremy L Thompson           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
2944b3e95d5SJeremy L Thompson 
2954b3e95d5SJeremy L Thompson           if (!has_backend_strides) {
2964b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
2974b3e95d5SJeremy L Thompson           }
2984b3e95d5SJeremy L Thompson           code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
299*f815fac9SJeremy L Thompson           code << "    ReadLVecStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << ","
3004b3e95d5SJeremy L Thompson                << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
301*f815fac9SJeremy L Thompson           break;
302*f815fac9SJeremy L Thompson         }
303*f815fac9SJeremy L Thompson         // LCOV_EXCL_START
304*f815fac9SJeremy L Thompson         case CEED_RESTRICTION_ORIENTED:
305*f815fac9SJeremy L Thompson         case CEED_RESTRICTION_CURL_ORIENTED:
306*f815fac9SJeremy L Thompson         case CEED_RESTRICTION_POINTS:
307*f815fac9SJeremy L Thompson           break;  // TODO: Not implemented
308*f815fac9SJeremy L Thompson                   // LCOV_EXCL_STOP
3094b3e95d5SJeremy L Thompson       }
3104b3e95d5SJeremy L Thompson     }
3114b3e95d5SJeremy L Thompson   } else {
3124b3e95d5SJeremy L Thompson     // Output
313*f815fac9SJeremy L Thompson     switch (rstr_type) {
314*f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STANDARD: {
3154b3e95d5SJeremy L Thompson         CeedInt comp_stride;
3164b3e95d5SJeremy L Thompson 
3174b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3184b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3194b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3204b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
3214b3e95d5SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
322*f815fac9SJeremy L Thompson         code << "    WriteLVecStandard" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size"
323*f815fac9SJeremy L Thompson              << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
324*f815fac9SJeremy L Thompson         break;
325*f815fac9SJeremy L Thompson       }
326*f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STRIDED: {
3274b3e95d5SJeremy L Thompson         bool    has_backend_strides;
3284b3e95d5SJeremy L Thompson         CeedInt num_elem;
3294b3e95d5SJeremy L Thompson 
3304b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
3314b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3324b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3334b3e95d5SJeremy L Thompson 
3344b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
3354b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
3364b3e95d5SJeremy L Thompson         }
3374b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
338*f815fac9SJeremy L Thompson         code << "    WriteLVecStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << ","
3394b3e95d5SJeremy L Thompson              << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
340*f815fac9SJeremy L Thompson         break;
341*f815fac9SJeremy L Thompson       }
342*f815fac9SJeremy L Thompson       // LCOV_EXCL_START
343*f815fac9SJeremy L Thompson       case CEED_RESTRICTION_ORIENTED:
344*f815fac9SJeremy L Thompson       case CEED_RESTRICTION_CURL_ORIENTED:
345*f815fac9SJeremy L Thompson       case CEED_RESTRICTION_POINTS:
346*f815fac9SJeremy L Thompson         break;  // TODO: Not implemented
347*f815fac9SJeremy L Thompson                 // LCOV_EXCL_STOP
3484b3e95d5SJeremy L Thompson     }
3494b3e95d5SJeremy L Thompson   }
3504b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
3514b3e95d5SJeremy L Thompson }
3524b3e95d5SJeremy L Thompson 
3534b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
3544b3e95d5SJeremy L Thompson // Basis
3554b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
3564b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
3574b3e95d5SJeremy L Thompson                                                 CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input,
3584b3e95d5SJeremy L Thompson                                                 bool use_3d_slices) {
3594b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
3604b3e95d5SJeremy L Thompson   std::string         P_name = "P_1d" + var_suffix, Q_name = "Q_1d";
3614b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
3624b3e95d5SJeremy L Thompson   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
3634b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
3644b3e95d5SJeremy L Thompson   CeedBasis           basis;
3654b3e95d5SJeremy L Thompson 
3664b3e95d5SJeremy L Thompson   // Get field data
3674b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
3684b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
3694b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
3704b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
3714b3e95d5SJeremy L Thompson   }
3724b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
3734b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
3744b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
3754b3e95d5SJeremy L Thompson   }
3764b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
3774b3e95d5SJeremy L Thompson 
3784b3e95d5SJeremy L Thompson   // Basis
3794b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
3804b3e95d5SJeremy L Thompson   if (is_input) {
3814b3e95d5SJeremy L Thompson     switch (eval_mode) {
3824b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
3834b3e95d5SJeremy L Thompson         if (!use_3d_slices) {
3844b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
3854b3e95d5SJeremy L Thompson         }
3864b3e95d5SJeremy L Thompson         break;
3874b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
3884b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
3894b3e95d5SJeremy L Thompson         code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
3904b3e95d5SJeremy L Thompson              << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
3914b3e95d5SJeremy L Thompson         break;
3924b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
3934b3e95d5SJeremy L Thompson         if (use_3d_slices) {
3944b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
3954b3e95d5SJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
3964b3e95d5SJeremy L Thompson                << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
3974b3e95d5SJeremy L Thompson         } else {
3984b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n";
3994b3e95d5SJeremy L Thompson           code << "    Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" << var_suffix
4004b3e95d5SJeremy L Thompson                << ", P_1d" << var_suffix << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q"
4014b3e95d5SJeremy L Thompson                << var_suffix << ");\n";
4024b3e95d5SJeremy L Thompson         }
4034b3e95d5SJeremy L Thompson         break;
4044b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
4054b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
4064b3e95d5SJeremy L Thompson 
4074b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[" << Q_name << "];\n";
4084b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
4094b3e95d5SJeremy L Thompson         data->W = basis_data->d_q_weight_1d;
4104b3e95d5SJeremy L Thompson         code << "    Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
4114b3e95d5SJeremy L Thompson         break;
4124b3e95d5SJeremy L Thompson       }
4134b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
4144b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
4154b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
4164b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
4174b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
4184b3e95d5SJeremy L Thompson     }
4194b3e95d5SJeremy L Thompson   } else {
4204b3e95d5SJeremy L Thompson     switch (eval_mode) {
4214b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
4224b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
4234b3e95d5SJeremy L Thompson         break;  // No action
4244b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
425e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
4264b3e95d5SJeremy L Thompson         code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
4274b3e95d5SJeremy L Thompson              << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
4284b3e95d5SJeremy L Thompson         break;
4294b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
430e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
4314b3e95d5SJeremy L Thompson         if (use_3d_slices) {
4324b3e95d5SJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
4334b3e95d5SJeremy L Thompson                << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
4344b3e95d5SJeremy L Thompson         } else {
4354b3e95d5SJeremy L Thompson           code << "    GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp"
4364b3e95d5SJeremy L Thompson                << var_suffix << ", " << P_name << "," << Q_name << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix
4374b3e95d5SJeremy L Thompson                << ", r_e" << var_suffix << ");\n";
4384b3e95d5SJeremy L Thompson         }
4394b3e95d5SJeremy L Thompson         break;
4404b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
4414b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
4424b3e95d5SJeremy L Thompson         break;  // Should not occur
4434b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
4444b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
4454b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
4464b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
4474b3e95d5SJeremy L Thompson     }
4484b3e95d5SJeremy L Thompson   }
4494b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4504b3e95d5SJeremy L Thompson }
4514b3e95d5SJeremy L Thompson 
4524b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4534b3e95d5SJeremy L Thompson // QFunction
4544b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4554b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt dim, CeedInt num_input_fields,
4564b3e95d5SJeremy L Thompson                                                     CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields,
4574b3e95d5SJeremy L Thompson                                                     CeedInt num_output_fields, CeedOperatorField *op_output_fields,
4584b3e95d5SJeremy L Thompson                                                     CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d,
4594b3e95d5SJeremy L Thompson                                                     bool use_3d_slices) {
4604b3e95d5SJeremy L Thompson   std::string         Q_name    = "Q_1d";
4614b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
4624b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
4634b3e95d5SJeremy L Thompson 
4644b3e95d5SJeremy L Thompson   // Setup output arays
4654b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
4664b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
4674b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
4684b3e95d5SJeremy L Thompson 
4694b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
4704b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
4714b3e95d5SJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) {
4724b3e95d5SJeremy L Thompson       code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
4734b3e95d5SJeremy L Thompson     }
4744b3e95d5SJeremy L Thompson     if (eval_mode == CEED_EVAL_GRAD) {
4754b3e95d5SJeremy L Thompson       if (use_3d_slices) {
4764b3e95d5SJeremy L Thompson         // Accumulator for gradient slices
4774b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
4784b3e95d5SJeremy L Thompson         code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
4794b3e95d5SJeremy L Thompson         code << "      r_q" << var_suffix << "[i] = 0.0;\n";
4804b3e95d5SJeremy L Thompson         code << "    }\n";
4814b3e95d5SJeremy L Thompson       } else {
4824b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n";
4834b3e95d5SJeremy L Thompson       }
4844b3e95d5SJeremy L Thompson     }
4854b3e95d5SJeremy L Thompson   }
4864b3e95d5SJeremy L Thompson 
4874b3e95d5SJeremy L Thompson   // We treat quadrature points per slice in 3d to save registers
4884b3e95d5SJeremy L Thompson   if (use_3d_slices) {
4894b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
4904b3e95d5SJeremy L Thompson     code << "    #pragma unroll\n";
4914b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
4924b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
4934b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
4944b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
4954b3e95d5SJeremy L Thompson 
4964b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
4974b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
4984b3e95d5SJeremy L Thompson       // Basis action
4994b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
5004b3e95d5SJeremy L Thompson       switch (eval_mode) {
5014b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
5024b3e95d5SJeremy L Thompson           bool is_strided;
5034b3e95d5SJeremy L Thompson 
5044b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
5054b3e95d5SJeremy L Thompson 
5064b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
5074b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
5084b3e95d5SJeremy L Thompson           if (is_strided) {
5094b3e95d5SJeremy L Thompson             bool    has_backend_strides;
5104b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
5114b3e95d5SJeremy L Thompson 
5124b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
5134b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
5144b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
5154b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
5164b3e95d5SJeremy L Thompson 
5174b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
5184b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
5194b3e95d5SJeremy L Thompson             }
5204b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
521*f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << "," << strides[0] << "," << strides[1] << ","
5224b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
5234b3e95d5SJeremy L Thompson           } else {
5244b3e95d5SJeremy L Thompson             CeedSize                 l_size = 0;
5254b3e95d5SJeremy L Thompson             CeedInt                  comp_stride;
5264b3e95d5SJeremy L Thompson             CeedElemRestriction_Hip *rstr_data;
5274b3e95d5SJeremy L Thompson 
5284b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
5294b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
5304b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
5314b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
5324b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
5334b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
534*f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
5354b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
5364b3e95d5SJeremy L Thompson           }
5374b3e95d5SJeremy L Thompson           break;
5384b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
5394b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
5404b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
5414b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
5424b3e95d5SJeremy L Thompson           code << "      }\n";
5434b3e95d5SJeremy L Thompson           break;
5444b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
5454b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
546*f815fac9SJeremy L Thompson           code << "      GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix
547*f815fac9SJeremy L Thompson                << ", r_s" << var_suffix << ");\n";
5484b3e95d5SJeremy L Thompson           break;
5494b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
5504b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
5514b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
5524b3e95d5SJeremy L Thompson           break;  // No action
5534b3e95d5SJeremy L Thompson                   // LCOV_EXCL_START
5544b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
5554b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
5564b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
5574b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
5584b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
5594b3e95d5SJeremy L Thompson       }
5604b3e95d5SJeremy L Thompson     }
5614b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
5624b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
5634b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
5644b3e95d5SJeremy L Thompson 
5654b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
5664b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
5674b3e95d5SJeremy L Thompson       // Basis action
5684b3e95d5SJeremy L Thompson       switch (eval_mode) {
5694b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
5704b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
5714b3e95d5SJeremy L Thompson           break;  // No action
5724b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
5734b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
5744b3e95d5SJeremy L Thompson           break;
5754b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
5764b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
5774b3e95d5SJeremy L Thompson           break;
5784b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
5794b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
5804b3e95d5SJeremy L Thompson           break;  // Should not occur
5814b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
5824b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
5834b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
5844b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
5854b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
5864b3e95d5SJeremy L Thompson       }
5874b3e95d5SJeremy L Thompson     }
5884b3e95d5SJeremy L Thompson   } else {
5894b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
5904b3e95d5SJeremy L Thompson     code << "    {\n";
5914b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
5924b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
5934b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
5944b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
5954b3e95d5SJeremy L Thompson     }
5964b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
5974b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
5984b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
5994b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
6004b3e95d5SJeremy L Thompson     }
6014b3e95d5SJeremy L Thompson   }
6024b3e95d5SJeremy L Thompson 
6034b3e95d5SJeremy L Thompson   // Input and output buffers
6044b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
6054b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
6064b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
6074b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
6084b3e95d5SJeremy L Thompson     code << "      // ------ Input field " << i << "\n";
6094b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
6104b3e95d5SJeremy L Thompson   }
6114b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
6124b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
6134b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
6144b3e95d5SJeremy L Thompson     code << "      // ------ Output field " << i << "\n";
6154b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
6164b3e95d5SJeremy L Thompson   }
6174b3e95d5SJeremy L Thompson 
6184b3e95d5SJeremy L Thompson   // Apply QFunction
6194b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
6204b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
6214b3e95d5SJeremy L Thompson   if (dim != 3 || use_3d_slices) {
6224b3e95d5SJeremy L Thompson     code << "1";
6234b3e95d5SJeremy L Thompson   } else {
6244b3e95d5SJeremy L Thompson     code << "Q_1d";
6254b3e95d5SJeremy L Thompson   }
6264b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
6274b3e95d5SJeremy L Thompson 
6284b3e95d5SJeremy L Thompson   // Copy or apply transpose grad, if needed
6294b3e95d5SJeremy L Thompson   if (use_3d_slices) {
6304b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
6314b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
6324b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
6334b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
6344b3e95d5SJeremy L Thompson 
6354b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
6364b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
6374b3e95d5SJeremy L Thompson       // Basis action
6384b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6394b3e95d5SJeremy L Thompson       switch (eval_mode) {
6404b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
6414b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
6424b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
6434b3e95d5SJeremy L Thompson           code << "      }\n";
6444b3e95d5SJeremy L Thompson           break;  // No action
6454b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
6464b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
6474b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
6484b3e95d5SJeremy L Thompson           code << "      }\n";
6494b3e95d5SJeremy L Thompson           break;
6504b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
651*f815fac9SJeremy L Thompson           code << "      GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G"
652*f815fac9SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
6534b3e95d5SJeremy L Thompson           break;
6544b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
6554b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
6564b3e95d5SJeremy L Thompson           break;  // Should not occur
6574b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
6584b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
6594b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
6604b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
6614b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
6624b3e95d5SJeremy L Thompson       }
6634b3e95d5SJeremy L Thompson     }
6644b3e95d5SJeremy L Thompson   }
6654b3e95d5SJeremy L Thompson   code << "    }\n";
6664b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
6674b3e95d5SJeremy L Thompson }
6684b3e95d5SJeremy L Thompson 
6694b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
6709e201c85SYohann // Build single operator kernel
6717d8d0e25Snbeams //------------------------------------------------------------------------------
672eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op) {
6734b3e95d5SJeremy L Thompson   bool                   is_tensor = true, use_3d_slices = false;
6747d8d0e25Snbeams   Ceed                   ceed;
6754b3e95d5SJeremy L Thompson   CeedInt                Q_1d, num_input_fields, num_output_fields, dim = 1;
676b7453713SJeremy L Thompson   CeedQFunctionField    *qf_input_fields, *qf_output_fields;
677b7453713SJeremy L Thompson   CeedQFunction_Hip_gen *qf_data;
678b7453713SJeremy L Thompson   CeedQFunction          qf;
679b7453713SJeremy L Thompson   CeedOperatorField     *op_input_fields, *op_output_fields;
680b7453713SJeremy L Thompson   CeedOperator_Hip_gen  *data;
6814b3e95d5SJeremy L Thompson   std::ostringstream     code;
6824b3e95d5SJeremy L Thompson 
6834b3e95d5SJeremy L Thompson   {
6844b3e95d5SJeremy L Thompson     bool is_setup_done;
685b7453713SJeremy L Thompson 
686b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
687b7453713SJeremy L Thompson     if (is_setup_done) return CEED_ERROR_SUCCESS;
6884b3e95d5SJeremy L Thompson   }
689b7453713SJeremy L Thompson 
690b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
691b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
692b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
693b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
694b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
695b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
6967d8d0e25Snbeams 
6974b3e95d5SJeremy L Thompson   // Get operator data
6984b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelData_Hip_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
6994b3e95d5SJeremy L Thompson                                                       qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
7004b3e95d5SJeremy L Thompson   if (dim == 0) dim = 1;
7014b3e95d5SJeremy L Thompson   data->dim = dim;
7024b3e95d5SJeremy L Thompson   if (Q_1d == 0) {
7034b3e95d5SJeremy L Thompson     CeedInt Q;
7044b3e95d5SJeremy L Thompson 
7054b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
7064b3e95d5SJeremy L Thompson     Q_1d = Q;
7074b3e95d5SJeremy L Thompson   }
7084b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
7094b3e95d5SJeremy L Thompson 
7100b454692Sjeremylt   // Check for restriction only identity operator
7114b3e95d5SJeremy L Thompson   {
7124b3e95d5SJeremy L Thompson     bool is_identity_qf;
7134b3e95d5SJeremy L Thompson 
7142b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
7150b454692Sjeremylt     if (is_identity_qf) {
7169e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
717b7453713SJeremy L Thompson 
7182b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
7192b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
7206574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
7216574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
7220b454692Sjeremylt     }
7234b3e95d5SJeremy L Thompson   }
724b2165e7aSSebastian Grimberg 
725b2165e7aSSebastian Grimberg   // Load basis source files
7264b3e95d5SJeremy L Thompson   // TODO: Add non-tensor, AtPoints
7279c25dd66SJeremy L Thompson   code << "// Tensor basis source\n";
7289c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-templates.h>\n\n";
7299c25dd66SJeremy L Thompson   code << "// CodeGen operator source\n";
7309c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/hip/hip-gen-templates.h>\n\n";
7317d8d0e25Snbeams 
7324b3e95d5SJeremy L Thompson   // Get QFunction name
7334b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
7344b3e95d5SJeremy L Thompson   std::string operator_name;
7354b3e95d5SJeremy L Thompson 
73609095acaSJeremy L Thompson   operator_name = "CeedKernelHipGenOperator_" + qfunction_name;
7377d8d0e25Snbeams 
7389e201c85SYohann   // Define CEED_Q_VLA
7399e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
7404b3e95d5SJeremy L Thompson   if (dim != 3 || use_3d_slices) {
7419e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
7429e201c85SYohann   } else {
7439e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
7449e201c85SYohann   }
7459e201c85SYohann 
7464b3e95d5SJeremy L Thompson   // Add user QFunction source
7474b3e95d5SJeremy L Thompson   {
7489c25dd66SJeremy L Thompson     const char *source_path;
7494b3e95d5SJeremy L Thompson 
7509c25dd66SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
7519c25dd66SJeremy L Thompson     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/hip/gen backend requires QFunction source code file");
7529c25dd66SJeremy L Thompson 
7539c25dd66SJeremy L Thompson     code << "// User QFunction source\n";
7549c25dd66SJeremy L Thompson     code << "#include \"" << source_path << "\"\n\n";
7554b3e95d5SJeremy L Thompson   }
7567d8d0e25Snbeams 
7577d8d0e25Snbeams   // Setup
7587d8d0e25Snbeams   code << "\n// -----------------------------------------------------------------------------\n";
7594b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
7604b3e95d5SJeremy L Thompson   code << "// \n";
7614b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
7624b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
7634b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
7644b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice  vector register\n";
7654b3e95d5SJeremy L Thompson   code << "// \n";
7664b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
7674b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
7684b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
769b3e1519bSnbeams   code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
7702b730f8bSJeremy L Thompson   code << "__global__ void " << operator_name
7712b730f8bSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Hip indices, Fields_Hip fields, Fields_Hip B, Fields_Hip G, CeedScalar* W) {\n";
7724b3e95d5SJeremy L Thompson 
7734b3e95d5SJeremy L Thompson   // Scratch buffers
7749e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
7754b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
7764b3e95d5SJeremy L Thompson 
7772b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
7789e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
7794b3e95d5SJeremy L Thompson       code << "  const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n";
7807d8d0e25Snbeams     }
7817d8d0e25Snbeams   }
7829e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
7834b3e95d5SJeremy L Thompson     code << "  CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n";
7847d8d0e25Snbeams   }
7857d8d0e25Snbeams 
7869e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
7879e201c85SYohann   code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
7887d8d0e25Snbeams 
7894b3e95d5SJeremy L Thompson   // Shared data
7904b3e95d5SJeremy L Thompson   code << "  extern __shared__ CeedScalar slice[];\n";
7919e201c85SYohann   code << "  SharedData_Hip data;\n";
7929e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
7939e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
7949e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
7959e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
7969e201c85SYohann   code << "  data.slice = slice + data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n";
7977d8d0e25Snbeams 
7987d8d0e25Snbeams   // Initialize constants, and matrices B and G
7994b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
8009e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
8015a5594ffSJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_input_fields[i], qf_input_fields[i], Q_1d, true, use_3d_slices));
8027d8d0e25Snbeams   }
8034b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
8049e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
8055a5594ffSJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices));
8064b3e95d5SJeremy L Thompson   }
8077d8d0e25Snbeams 
8084b3e95d5SJeremy L Thompson   // Loop over all elements
8094b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
8107d8d0e25Snbeams   code << "  __syncthreads();\n";
8119e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
8124b3e95d5SJeremy L Thompson 
813e93651e5SJeremy L Thompson   // -- Compute minimum buffer space needed
814e93651e5SJeremy L Thompson   CeedInt max_rstr_buffer_size = 0;
815e93651e5SJeremy L Thompson 
816e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
817e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
818e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
819e93651e5SJeremy L Thompson 
820e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
821e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
822e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
823e93651e5SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size);
824681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
825e93651e5SJeremy L Thompson   }
826e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
827e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
828e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
829e93651e5SJeremy L Thompson 
830e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
831e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
832e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
833e93651e5SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size);
834681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
835e93651e5SJeremy L Thompson   }
836e93651e5SJeremy L Thompson   code << "    // Scratch restriction buffer space\n";
837e93651e5SJeremy L Thompson   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
838e93651e5SJeremy L Thompson 
839e93651e5SJeremy L Thompson   // -- Determine best input field processing order
840e93651e5SJeremy L Thompson   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
841e93651e5SJeremy L Thompson 
842e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
843e93651e5SJeremy L Thompson     field_rstr_in_buffer[i] = -1;
844e93651e5SJeremy L Thompson     input_field_order[i]    = -1;
845e93651e5SJeremy L Thompson   }
846e93651e5SJeremy L Thompson   {
847e93651e5SJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
848e93651e5SJeremy L Thompson     CeedInt curr_index = 0;
849e93651e5SJeremy L Thompson 
850e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
851e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
852e93651e5SJeremy L Thompson       CeedVector          vec_i;
853e93651e5SJeremy L Thompson       CeedElemRestriction rstr_i;
854e93651e5SJeremy L Thompson 
855e93651e5SJeremy L Thompson       if (is_ordered[i]) continue;
856e93651e5SJeremy L Thompson       field_rstr_in_buffer[i]       = i;
857e93651e5SJeremy L Thompson       is_ordered[i]                 = true;
858e93651e5SJeremy L Thompson       input_field_order[curr_index] = i;
859e93651e5SJeremy L Thompson       curr_index++;
860034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
861e93651e5SJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
862e93651e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
863e93651e5SJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
864e93651e5SJeremy L Thompson         CeedVector          vec_j;
865e93651e5SJeremy L Thompson         CeedElemRestriction rstr_j;
866e93651e5SJeremy L Thompson 
867e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
868e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
869e93651e5SJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
870e93651e5SJeremy L Thompson           field_rstr_in_buffer[j]       = i;
871e93651e5SJeremy L Thompson           is_ordered[j]                 = true;
872e93651e5SJeremy L Thompson           input_field_order[curr_index] = j;
873e93651e5SJeremy L Thompson           curr_index++;
874e93651e5SJeremy L Thompson         }
875e93651e5SJeremy L Thompson       }
876e93651e5SJeremy L Thompson     }
877e93651e5SJeremy L Thompson   }
878e93651e5SJeremy L Thompson 
8794b3e95d5SJeremy L Thompson   // -- Input restriction and basis
8804b3e95d5SJeremy L Thompson   code << "    // -- Input field restrictions and basis actions\n";
8819e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
882e93651e5SJeremy L Thompson     CeedInt f = input_field_order[i];
883e93651e5SJeremy L Thompson 
884e93651e5SJeremy L Thompson     code << "    // ---- Input field " << f << "\n";
8857d8d0e25Snbeams 
8864b3e95d5SJeremy L Thompson     // ---- Restriction
887e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], Q_1d,
888e93651e5SJeremy L Thompson                                                                true, use_3d_slices));
889b7453713SJeremy L Thompson 
8904b3e95d5SJeremy L Thompson     // ---- Basis action
891e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, use_3d_slices));
8927d8d0e25Snbeams   }
8937d8d0e25Snbeams 
8944b3e95d5SJeremy L Thompson   // -- Q function
8954b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Hip_gen(code, data, dim, num_input_fields, op_input_fields, qf_input_fields, num_output_fields,
8964b3e95d5SJeremy L Thompson                                                            op_output_fields, qf_output_fields, qfunction_name, Q_1d, use_3d_slices));
8977d8d0e25Snbeams 
8984b3e95d5SJeremy L Thompson   // -- Output basis and restriction
8994b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
9009e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
9014b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
902b7453713SJeremy L Thompson 
9034b3e95d5SJeremy L Thompson     // ---- Basis action
9044b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices));
9057d8d0e25Snbeams 
9064b3e95d5SJeremy L Thompson     // ---- Restriction
9074b3e95d5SJeremy L Thompson     CeedCallBackend(
908e93651e5SJeremy L Thompson         CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices));
9097d8d0e25Snbeams   }
9107d8d0e25Snbeams 
9114b3e95d5SJeremy L Thompson   // Close loop and function
9127d8d0e25Snbeams   code << "  }\n";
9137d8d0e25Snbeams   code << "}\n";
9147d8d0e25Snbeams   code << "// -----------------------------------------------------------------------------\n\n";
9157d8d0e25Snbeams 
916539ec17dSJeremy L Thompson   CeedInt block_sizes[3] = {0, 0, 0};
9179e201c85SYohann   CeedInt num_elem;
918b7453713SJeremy L Thompson 
9192b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
9202b730f8bSJeremy L Thompson   CeedCallBackend(BlockGridCalculate_Hip_gen(dim, num_elem, data->max_P_1d, Q_1d, block_sizes));
921eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE",
9222b730f8bSJeremy L Thompson                                   block_sizes[0] * block_sizes[1] * block_sizes[2]));
923eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op));
9242b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
9259bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
926c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
927e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
9287d8d0e25Snbeams }
9292a86cc9dSSebastian Grimberg 
9307d8d0e25Snbeams //------------------------------------------------------------------------------
931