xref: /libCEED/rust/libceed-sys/c-src/backends/hip-gen/ceed-hip-gen-operator-build.cpp (revision 3a2968d63a7f2ece086fcd3a62875aca8b9498aa)
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) {
28*3a2968d6SJeremy L Thompson   const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
29b3e1519bSnbeams   if (dim == 1) {
30*3a2968d6SJeremy L Thompson     CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
31b7453713SJeremy L Thompson 
329e201c85SYohann     elems_per_block = elems_per_block > 0 ? elems_per_block : 1;
33*3a2968d6SJeremy L Thompson     block_sizes[0]  = thread_1d;
34b3e1519bSnbeams     block_sizes[1]  = 1;
359e201c85SYohann     block_sizes[2]  = elems_per_block;
36b3e1519bSnbeams   } else if (dim == 2) {
37*3a2968d6SJeremy L Thompson     const CeedInt elems_per_block = thread_1d < 4 ? 16 : 2;
38b7453713SJeremy L Thompson 
39*3a2968d6SJeremy L Thompson     block_sizes[0] = thread_1d;
40*3a2968d6SJeremy L Thompson     block_sizes[1] = thread_1d;
419e201c85SYohann     block_sizes[2] = elems_per_block;
42b3e1519bSnbeams   } else if (dim == 3) {
43*3a2968d6SJeremy L Thompson     const CeedInt elems_per_block = thread_1d < 6 ? 4 : (thread_1d < 8 ? 2 : 1);
44b7453713SJeremy L Thompson 
45*3a2968d6SJeremy L Thompson     block_sizes[0] = thread_1d;
46*3a2968d6SJeremy L Thompson     block_sizes[1] = thread_1d;
479e201c85SYohann     block_sizes[2] = elems_per_block;
48b3e1519bSnbeams   }
49b3e1519bSnbeams   return CEED_ERROR_SUCCESS;
50b3e1519bSnbeams }
51b3e1519bSnbeams 
527d8d0e25Snbeams //------------------------------------------------------------------------------
534b3e95d5SJeremy L Thompson // Determine type of operator
544b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
554b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Hip_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
564b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
574b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor,
584b3e95d5SJeremy L Thompson                                                bool *use_3d_slices) {
594b3e95d5SJeremy L Thompson   // Find dim, P_1d, Q_1d
604b3e95d5SJeremy L Thompson   *max_P_1d  = 0;
614b3e95d5SJeremy L Thompson   *Q_1d      = 0;
624b3e95d5SJeremy L Thompson   *dim       = 0;
634b3e95d5SJeremy L Thompson   *is_tensor = true;
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     }
85*3a2968d6SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
864b3e95d5SJeremy L Thompson   }
874b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
884b3e95d5SJeremy L Thompson     CeedBasis basis;
894b3e95d5SJeremy L Thompson 
904b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
914b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
924b3e95d5SJeremy L Thompson       bool    is_field_tensor;
934b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
944b3e95d5SJeremy L Thompson 
954b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
964b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
974b3e95d5SJeremy L Thompson       CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
984b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
994b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
1004b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
1014b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
1024b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
1034b3e95d5SJeremy L Thompson       *dim = field_dim;
1044b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
1054b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
1064b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
1074b3e95d5SJeremy L Thompson     }
108*3a2968d6SJeremy L Thompson     CeedCallBackend(CeedBasisDestroy(&basis));
1094b3e95d5SJeremy L Thompson   }
1104b3e95d5SJeremy L Thompson 
1114b3e95d5SJeremy L Thompson   // Only use 3D collocated gradient parallelization strategy when gradient is computed
1124b3e95d5SJeremy L Thompson   *use_3d_slices = false;
1134b3e95d5SJeremy L Thompson   if (*dim == 3) {
1144b3e95d5SJeremy L Thompson     bool was_grad_found = false;
1154b3e95d5SJeremy L Thompson 
1164b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1174b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1184b3e95d5SJeremy L Thompson 
1194b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1204b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1214b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
1224b3e95d5SJeremy L Thompson         CeedBasis             basis;
1234b3e95d5SJeremy L Thompson 
1244b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
1254b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1264b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1274b3e95d5SJeremy L Thompson         was_grad_found = true;
128*3a2968d6SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1294b3e95d5SJeremy L Thompson       }
1304b3e95d5SJeremy L Thompson     }
1314b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
1324b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
1334b3e95d5SJeremy L Thompson 
1344b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1354b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
1364b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
1374b3e95d5SJeremy L Thompson         CeedBasis             basis;
1384b3e95d5SJeremy L Thompson 
1394b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
1404b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1414b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
1424b3e95d5SJeremy L Thompson         was_grad_found = true;
143*3a2968d6SJeremy L Thompson         CeedCallBackend(CeedBasisDestroy(&basis));
1444b3e95d5SJeremy L Thompson       }
1454b3e95d5SJeremy L Thompson     }
1464b3e95d5SJeremy L Thompson   }
1474b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1484b3e95d5SJeremy L Thompson }
1494b3e95d5SJeremy L Thompson 
1504b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1514b3e95d5SJeremy L Thompson // Setup fields
1524b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
1534b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedOperatorField op_field,
154*3a2968d6SJeremy L Thompson                                                     CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_at_points, bool use_3d_slices) {
1554b3e95d5SJeremy L Thompson   std::string           var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
1564b3e95d5SJeremy L Thompson   std::string           P_name = "P_1d" + var_suffix, Q_name = "Q_1d";
1574b3e95d5SJeremy L Thompson   std::string           option_name = (is_input ? "inputs" : "outputs");
1584b3e95d5SJeremy L Thompson   CeedEvalMode          eval_mode   = CEED_EVAL_NONE;
1594b3e95d5SJeremy L Thompson   CeedInt               elem_size = 0, num_comp = 0, P_1d = 0;
1604b3e95d5SJeremy L Thompson   CeedElemRestriction   elem_rstr;
1614b3e95d5SJeremy L Thompson   CeedBasis_Hip_shared *basis_data;
1624b3e95d5SJeremy L Thompson   CeedBasis             basis;
1634b3e95d5SJeremy L Thompson 
1644b3e95d5SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n";
1654b3e95d5SJeremy L Thompson 
1664b3e95d5SJeremy L Thompson   // Get field data
1674b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
1684b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
1694b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1704b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1714b3e95d5SJeremy L Thompson   }
172*3a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1734b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
1744b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
1754b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
1764b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
1774b3e95d5SJeremy L Thompson   }
1784b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
1794b3e95d5SJeremy L Thompson 
1804b3e95d5SJeremy L Thompson   // Set field constants
1814b3e95d5SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
1824b3e95d5SJeremy L Thompson     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
1834b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
1844b3e95d5SJeremy L Thompson   }
1854b3e95d5SJeremy L Thompson 
1864b3e95d5SJeremy L Thompson   // Load basis data
1874b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
1884b3e95d5SJeremy L Thompson   switch (eval_mode) {
1894b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
1904b3e95d5SJeremy L Thompson       break;
1914b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
192*3a2968d6SJeremy L Thompson       if (is_at_points) {
193*3a2968d6SJeremy L Thompson         // AtPoints
194*3a2968d6SJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
195*3a2968d6SJeremy L Thompson           CeedSize    interp_bytes;
196*3a2968d6SJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
197*3a2968d6SJeremy L Thompson 
198*3a2968d6SJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
199*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
200*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
201*3a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
202*3a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis),
203*3a2968d6SJeremy L Thompson                       hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice));
204*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
205*3a2968d6SJeremy L Thompson         }
206*3a2968d6SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
207*3a2968d6SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
208*3a2968d6SJeremy L Thompson       } else {
209*3a2968d6SJeremy L Thompson         // Standard quadrature
2104b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2114b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
212*3a2968d6SJeremy L Thompson       }
2134b3e95d5SJeremy L Thompson       code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n";
214f815fac9SJeremy L Thompson       code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
2154b3e95d5SJeremy L Thompson       break;
2164b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
217*3a2968d6SJeremy L Thompson       if (is_at_points) {
218*3a2968d6SJeremy L Thompson         // AtPoints
219*3a2968d6SJeremy L Thompson         if (!basis_data->d_chebyshev_interp_1d) {
220*3a2968d6SJeremy L Thompson           CeedSize    interp_bytes;
221*3a2968d6SJeremy L Thompson           CeedScalar *chebyshev_interp_1d;
222*3a2968d6SJeremy L Thompson 
223*3a2968d6SJeremy L Thompson           interp_bytes = P_1d * Q_1d * sizeof(CeedScalar);
224*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d));
225*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d));
226*3a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes));
227*3a2968d6SJeremy L Thompson           CeedCallHip(CeedBasisReturnCeed(basis),
228*3a2968d6SJeremy L Thompson                       hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice));
229*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedFree(&chebyshev_interp_1d));
230*3a2968d6SJeremy L Thompson         }
231*3a2968d6SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d;
232*3a2968d6SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d;
233*3a2968d6SJeremy L Thompson       } else {
234*3a2968d6SJeremy L Thompson         // Standard quadrature
2354b3e95d5SJeremy L Thompson         if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
2364b3e95d5SJeremy L Thompson         else data->B.outputs[i] = basis_data->d_interp_1d;
237*3a2968d6SJeremy L Thompson       }
2384b3e95d5SJeremy L Thompson       code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n";
239f815fac9SJeremy L Thompson       code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
240*3a2968d6SJeremy L Thompson       if (is_at_points) break;  // No G mat for AtPoints
2414b3e95d5SJeremy L Thompson       if (use_3d_slices) {
2424b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
2434b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
2444b3e95d5SJeremy L Thompson         code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n";
245f815fac9SJeremy L Thompson         code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2464b3e95d5SJeremy L Thompson       } else {
2474b3e95d5SJeremy L Thompson         bool has_collo_grad = basis_data->d_collo_grad_1d;
2484b3e95d5SJeremy L Thompson 
2494b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2504b3e95d5SJeremy L Thompson         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
2514b3e95d5SJeremy L Thompson         if (has_collo_grad) {
2524b3e95d5SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n";
253f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2544b3e95d5SJeremy L Thompson         } else {
2554b3e95d5SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * P_1d << "];\n";
256f815fac9SJeremy L Thompson           code << "  LoadMatrix<" << P_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
2574b3e95d5SJeremy L Thompson         }
2584b3e95d5SJeremy L Thompson       }
2594b3e95d5SJeremy L Thompson       break;
2604b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
2614b3e95d5SJeremy L Thompson       break;  // No action
2624b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
2634b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
2644b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
2654b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
2664b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
2674b3e95d5SJeremy L Thompson   }
268*3a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
2694b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
2704b3e95d5SJeremy L Thompson }
2714b3e95d5SJeremy L Thompson 
2724b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2734b3e95d5SJeremy L Thompson // Restriction
2744b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
2754b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
276e93651e5SJeremy L Thompson                                                       CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field,
277*3a2968d6SJeremy L Thompson                                                       CeedInt Q_1d, bool is_input, bool is_at_points, bool use_3d_slices) {
2784b3e95d5SJeremy L Thompson   std::string              var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
2794b3e95d5SJeremy L Thompson   std::string              P_name     = "P_1d" + var_suffix;
2804b3e95d5SJeremy L Thompson   CeedEvalMode             eval_mode  = CEED_EVAL_NONE;
2814b3e95d5SJeremy L Thompson   CeedInt                  elem_size = 0, num_comp = 0, P_1d = 0;
2824b3e95d5SJeremy L Thompson   CeedSize                 l_size;
283f815fac9SJeremy L Thompson   CeedRestrictionType      rstr_type = CEED_RESTRICTION_STANDARD;
2844b3e95d5SJeremy L Thompson   CeedElemRestriction_Hip *rstr_data;
2854b3e95d5SJeremy L Thompson   CeedElemRestriction      elem_rstr;
2864b3e95d5SJeremy L Thompson   CeedBasis                basis;
2874b3e95d5SJeremy L Thompson 
2884b3e95d5SJeremy L Thompson   // Get field data
2894b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
2904b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
291f815fac9SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type));
2924b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
2934b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
2944b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
2954b3e95d5SJeremy L Thompson   }
2964b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
2974b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
2984b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
2994b3e95d5SJeremy L Thompson   }
300*3a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
3014b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
3024b3e95d5SJeremy L Thompson 
3034b3e95d5SJeremy L Thompson   // Restriction
3044b3e95d5SJeremy L Thompson   if (is_input) {
3054b3e95d5SJeremy L Thompson     // Input
306e93651e5SJeremy L Thompson     if (field_input_buffer[i] != i) {
307e93651e5SJeremy L Thompson       std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]);
308e93651e5SJeremy L Thompson 
309e93651e5SJeremy L Thompson       // Restriction was already done for previous input
310e93651e5SJeremy L Thompson       code << "    CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n";
311*3a2968d6SJeremy L Thompson     } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) {
312*3a2968d6SJeremy L Thompson       if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) {
313e93651e5SJeremy L Thompson         // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated
3144b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
315*3a2968d6SJeremy L Thompson       } else if (rstr_type != CEED_RESTRICTION_POINTS) {
316e93651e5SJeremy L Thompson         // Otherwise we're using the scratch space
317e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
318e93651e5SJeremy L Thompson       }
319f815fac9SJeremy L Thompson       switch (rstr_type) {
320f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STANDARD: {
3214b3e95d5SJeremy L Thompson           CeedInt comp_stride;
3224b3e95d5SJeremy L Thompson 
3234b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3244b3e95d5SJeremy L Thompson           code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3254b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3264b3e95d5SJeremy L Thompson           code << "    // CompStride: " << comp_stride << "\n";
3274b3e95d5SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
328f815fac9SJeremy L Thompson           code << "    ReadLVecStandard" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size"
329f815fac9SJeremy L Thompson                << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
330f815fac9SJeremy L Thompson           break;
331f815fac9SJeremy L Thompson         }
332f815fac9SJeremy L Thompson         case CEED_RESTRICTION_STRIDED: {
3334b3e95d5SJeremy L Thompson           bool    has_backend_strides;
3344b3e95d5SJeremy L Thompson           CeedInt num_elem;
3354b3e95d5SJeremy L Thompson 
3364b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
3374b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3384b3e95d5SJeremy L Thompson           CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3394b3e95d5SJeremy L Thompson 
3404b3e95d5SJeremy L Thompson           if (!has_backend_strides) {
3414b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
3424b3e95d5SJeremy L Thompson           }
3434b3e95d5SJeremy L Thompson           code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
344f815fac9SJeremy L Thompson           code << "    ReadLVecStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << ","
3454b3e95d5SJeremy L Thompson                << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
346f815fac9SJeremy L Thompson           break;
347f815fac9SJeremy L Thompson         }
348*3a2968d6SJeremy L Thompson         case CEED_RESTRICTION_POINTS: {
349*3a2968d6SJeremy L Thompson           CeedInt comp_stride;
350*3a2968d6SJeremy L Thompson 
351*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
352*3a2968d6SJeremy L Thompson           code << "    const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
353*3a2968d6SJeremy L Thompson           data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
354*3a2968d6SJeremy L Thompson           break;
355*3a2968d6SJeremy L Thompson         }
356f815fac9SJeremy L Thompson         // LCOV_EXCL_START
357f815fac9SJeremy L Thompson         case CEED_RESTRICTION_ORIENTED:
358f815fac9SJeremy L Thompson         case CEED_RESTRICTION_CURL_ORIENTED:
359f815fac9SJeremy L Thompson           break;  // TODO: Not implemented
360f815fac9SJeremy L Thompson                   // LCOV_EXCL_STOP
3614b3e95d5SJeremy L Thompson       }
3624b3e95d5SJeremy L Thompson     }
3634b3e95d5SJeremy L Thompson   } else {
3644b3e95d5SJeremy L Thompson     // Output
365f815fac9SJeremy L Thompson     switch (rstr_type) {
366f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STANDARD: {
3674b3e95d5SJeremy L Thompson         CeedInt comp_stride;
3684b3e95d5SJeremy L Thompson 
3694b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
3704b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
3714b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
3724b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
3734b3e95d5SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
374f815fac9SJeremy L Thompson         code << "    WriteLVecStandard" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size"
375f815fac9SJeremy L Thompson              << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
376f815fac9SJeremy L Thompson         break;
377f815fac9SJeremy L Thompson       }
378f815fac9SJeremy L Thompson       case CEED_RESTRICTION_STRIDED: {
3794b3e95d5SJeremy L Thompson         bool    has_backend_strides;
3804b3e95d5SJeremy L Thompson         CeedInt num_elem;
3814b3e95d5SJeremy L Thompson 
3824b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
3834b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
3844b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
3854b3e95d5SJeremy L Thompson 
3864b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
3874b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
3884b3e95d5SJeremy L Thompson         }
3894b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
390f815fac9SJeremy L Thompson         code << "    WriteLVecStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << ","
3914b3e95d5SJeremy L Thompson              << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
392f815fac9SJeremy L Thompson         break;
393f815fac9SJeremy L Thompson       }
394*3a2968d6SJeremy L Thompson       case CEED_RESTRICTION_POINTS:
395*3a2968d6SJeremy L Thompson         data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
396*3a2968d6SJeremy L Thompson         break;
397f815fac9SJeremy L Thompson       // LCOV_EXCL_START
398f815fac9SJeremy L Thompson       case CEED_RESTRICTION_ORIENTED:
399f815fac9SJeremy L Thompson       case CEED_RESTRICTION_CURL_ORIENTED:
400f815fac9SJeremy L Thompson         break;  // TODO: Not implemented
401f815fac9SJeremy L Thompson                 // LCOV_EXCL_STOP
4024b3e95d5SJeremy L Thompson     }
4034b3e95d5SJeremy L Thompson   }
404*3a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4054b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
4064b3e95d5SJeremy L Thompson }
4074b3e95d5SJeremy L Thompson 
4084b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4094b3e95d5SJeremy L Thompson // Basis
4104b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
4114b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
4124b3e95d5SJeremy L Thompson                                                 CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input,
413*3a2968d6SJeremy L Thompson                                                 bool is_at_points, bool use_3d_slices) {
4144b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
4154b3e95d5SJeremy L Thompson   std::string         P_name = "P_1d" + var_suffix, Q_name = "Q_1d";
4164b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
4174b3e95d5SJeremy L Thompson   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
4184b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
4194b3e95d5SJeremy L Thompson   CeedBasis           basis;
4204b3e95d5SJeremy L Thompson 
4214b3e95d5SJeremy L Thompson   // Get field data
4224b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
4234b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
4244b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
4254b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
4264b3e95d5SJeremy L Thompson   }
427*3a2968d6SJeremy L Thompson   CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
4284b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
4294b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
4304b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
4314b3e95d5SJeremy L Thompson   }
4324b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
4334b3e95d5SJeremy L Thompson 
4344b3e95d5SJeremy L Thompson   // Basis
4354b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
4364b3e95d5SJeremy L Thompson   if (is_input) {
4374b3e95d5SJeremy L Thompson     switch (eval_mode) {
4384b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
439*3a2968d6SJeremy L Thompson         if (!use_3d_slices && !is_at_points) {
4404b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
4414b3e95d5SJeremy L Thompson         }
4424b3e95d5SJeremy L Thompson         break;
4434b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
444*3a2968d6SJeremy L Thompson         if (is_at_points) {
445*3a2968d6SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim > 2 ? Q_name : "1") << "];\n";
446*3a2968d6SJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
447*3a2968d6SJeremy L Thompson                << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n";
448*3a2968d6SJeremy L Thompson         } else {
4494b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
4504b3e95d5SJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
4514b3e95d5SJeremy L Thompson                << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
452*3a2968d6SJeremy L Thompson         }
4534b3e95d5SJeremy L Thompson         break;
4544b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
455*3a2968d6SJeremy L Thompson         if (is_at_points) {
456*3a2968d6SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim > 2 ? Q_name : "1") << "];\n";
457*3a2968d6SJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
458*3a2968d6SJeremy L Thompson                << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n";
459*3a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
4604b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
4614b3e95d5SJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
4624b3e95d5SJeremy L Thompson                << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
4634b3e95d5SJeremy L Thompson         } else {
4644b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n";
4654b3e95d5SJeremy L Thompson           code << "    Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" << var_suffix
4664b3e95d5SJeremy L Thompson                << ", P_1d" << var_suffix << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q"
4674b3e95d5SJeremy L Thompson                << var_suffix << ");\n";
4684b3e95d5SJeremy L Thompson         }
4694b3e95d5SJeremy L Thompson         break;
4704b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
471*3a2968d6SJeremy L Thompson         if (is_at_points) {
472*3a2968d6SJeremy L Thompson           code << "    // Nothing to do AtPoints\n";
473*3a2968d6SJeremy L Thompson         } else {
4744b3e95d5SJeremy L Thompson           CeedBasis_Hip_shared *basis_data;
4754b3e95d5SJeremy L Thompson 
4764b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[" << Q_name << "];\n";
4774b3e95d5SJeremy L Thompson           CeedCallBackend(CeedBasisGetData(basis, &basis_data));
4784b3e95d5SJeremy L Thompson           data->W = basis_data->d_q_weight_1d;
4794b3e95d5SJeremy L Thompson           code << "    Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
480*3a2968d6SJeremy L Thompson         }
4814b3e95d5SJeremy L Thompson         break;
4824b3e95d5SJeremy L Thompson       }
4834b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
4844b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
4854b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
4864b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
4874b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
4884b3e95d5SJeremy L Thompson     }
4894b3e95d5SJeremy L Thompson   } else {
4904b3e95d5SJeremy L Thompson     switch (eval_mode) {
4914b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
4924b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
4934b3e95d5SJeremy L Thompson         break;  // No action
4944b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
495e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
496*3a2968d6SJeremy L Thompson         if (is_at_points) {
497*3a2968d6SJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
498*3a2968d6SJeremy L Thompson                << ">(data, r_c" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
499*3a2968d6SJeremy L Thompson         } else {
5004b3e95d5SJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
5014b3e95d5SJeremy L Thompson                << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
502*3a2968d6SJeremy L Thompson         }
5034b3e95d5SJeremy L Thompson         break;
5044b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
505e93651e5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n";
506*3a2968d6SJeremy L Thompson         if (is_at_points) {
507*3a2968d6SJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
508*3a2968d6SJeremy L Thompson                << ">(data, r_c" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
509*3a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
5104b3e95d5SJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
5114b3e95d5SJeremy L Thompson                << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
5124b3e95d5SJeremy L Thompson         } else {
5134b3e95d5SJeremy L Thompson           code << "    GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp"
5144b3e95d5SJeremy L Thompson                << var_suffix << ", " << P_name << "," << Q_name << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix
5154b3e95d5SJeremy L Thompson                << ", r_e" << var_suffix << ");\n";
5164b3e95d5SJeremy L Thompson         }
5174b3e95d5SJeremy L Thompson         break;
5184b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
5194b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
5204b3e95d5SJeremy L Thompson         break;  // Should not occur
5214b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
5224b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
5234b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
5244b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
5254b3e95d5SJeremy L Thompson     }
5264b3e95d5SJeremy L Thompson   }
527*3a2968d6SJeremy L Thompson   CeedCallBackend(CeedBasisDestroy(&basis));
5284b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
5294b3e95d5SJeremy L Thompson }
5304b3e95d5SJeremy L Thompson 
5314b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
5324b3e95d5SJeremy L Thompson // QFunction
5334b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
534*3a2968d6SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt dim, CeedInt max_num_points,
535*3a2968d6SJeremy L Thompson                                                     CeedInt num_input_fields, CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields,
5364b3e95d5SJeremy L Thompson                                                     CeedInt num_output_fields, CeedOperatorField *op_output_fields,
537*3a2968d6SJeremy L Thompson                                                     CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d, bool is_at_points,
5384b3e95d5SJeremy L Thompson                                                     bool use_3d_slices) {
5394b3e95d5SJeremy L Thompson   std::string         Q_name    = "Q_1d";
5404b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
5414b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
5424b3e95d5SJeremy L Thompson 
5438b97b69aSJeremy L Thompson   // Setup output arrays
5444b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
5454b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
5464b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
5474b3e95d5SJeremy L Thompson 
5484b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
5494b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
550*3a2968d6SJeremy L Thompson     switch (eval_mode) {
551*3a2968d6SJeremy L Thompson       case CEED_EVAL_NONE:
552*3a2968d6SJeremy L Thompson         if (is_at_points) {
553*3a2968d6SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n";
554*3a2968d6SJeremy L Thompson         } else {
5554b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
5564b3e95d5SJeremy L Thompson         }
557*3a2968d6SJeremy L Thompson         break;
558*3a2968d6SJeremy L Thompson       case CEED_EVAL_INTERP:
559*3a2968d6SJeremy L Thompson         if (is_at_points) {
560*3a2968d6SJeremy L Thompson           // Accumulator for point data
561*3a2968d6SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim > 2 ? Q_name : "1") << "];\n";
562*3a2968d6SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim > 2 ? Q_name : "1") << "; i++) {\n";
563*3a2968d6SJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
564*3a2968d6SJeremy L Thompson           code << "    }\n";
565*3a2968d6SJeremy L Thompson         } else {
566*3a2968d6SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
567*3a2968d6SJeremy L Thompson         }
568*3a2968d6SJeremy L Thompson         break;
569*3a2968d6SJeremy L Thompson       case CEED_EVAL_GRAD:
570*3a2968d6SJeremy L Thompson         if (is_at_points) {
571*3a2968d6SJeremy L Thompson           // Accumulator for point data
572*3a2968d6SJeremy L Thompson           code << "    CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim > 2 ? Q_name : "1") << "*dim];\n";
573*3a2968d6SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim > 2 ? Q_name : "1") << "; i++) {\n";
574*3a2968d6SJeremy L Thompson           code << "      r_c" << var_suffix << "[i] = 0.0;\n";
575*3a2968d6SJeremy L Thompson           code << "    }\n";
576*3a2968d6SJeremy L Thompson         } else if (use_3d_slices) {
5774b3e95d5SJeremy L Thompson           // Accumulator for gradient slices
5784b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
5794b3e95d5SJeremy L Thompson           code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
5804b3e95d5SJeremy L Thompson           code << "      r_q" << var_suffix << "[i] = 0.0;\n";
5814b3e95d5SJeremy L Thompson           code << "    }\n";
5824b3e95d5SJeremy L Thompson         } else {
5834b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n";
5844b3e95d5SJeremy L Thompson         }
585*3a2968d6SJeremy L Thompson         break;
586*3a2968d6SJeremy L Thompson       case CEED_EVAL_WEIGHT:
587*3a2968d6SJeremy L Thompson         break;
588*3a2968d6SJeremy L Thompson         // LCOV_EXCL_START
589*3a2968d6SJeremy L Thompson       case CEED_EVAL_DIV:
590*3a2968d6SJeremy L Thompson       case CEED_EVAL_CURL:
591*3a2968d6SJeremy L Thompson         break;  // TODO: Not implemented
592*3a2968d6SJeremy L Thompson                 // LCOV_EXCL_STOP
5934b3e95d5SJeremy L Thompson     }
5944b3e95d5SJeremy L Thompson   }
5954b3e95d5SJeremy L Thompson 
596*3a2968d6SJeremy L Thompson   if (is_at_points) {
597*3a2968d6SJeremy L Thompson     // We need to handle batches of points
598*3a2968d6SJeremy L Thompson     code << "\n    // Note: Using batches of points\n";
599*3a2968d6SJeremy L Thompson     code << "    const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n";
600*3a2968d6SJeremy L Thompson     code << "    #pragma unroll\n";
601*3a2968d6SJeremy L Thompson     code << "    for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n";
602*3a2968d6SJeremy L Thompson     code << "      const CeedInt p = i % max_num_points;\n\n";
603*3a2968d6SJeremy L Thompson 
604*3a2968d6SJeremy L Thompson     code << "      // -- Coordinates\n";
605*3a2968d6SJeremy L Thompson     code << "      CeedScalar r_x[dim];\n";
606*3a2968d6SJeremy 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";
607*3a2968d6SJeremy L Thompson 
608*3a2968d6SJeremy L Thompson     code << "      // -- Input fields\n";
609*3a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
610*3a2968d6SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
611*3a2968d6SJeremy L Thompson 
612*3a2968d6SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
613*3a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
614*3a2968d6SJeremy L Thompson       // Basis action
615*3a2968d6SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
616*3a2968d6SJeremy L Thompson       switch (eval_mode) {
617*3a2968d6SJeremy L Thompson         case CEED_EVAL_NONE:
618*3a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
619*3a2968d6SJeremy L Thompson           code << "      ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
620*3a2968d6SJeremy L Thompson                << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
621*3a2968d6SJeremy L Thompson           break;
622*3a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
623*3a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
624*3a2968d6SJeremy L Thompson           code << "      InterpAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
625*3a2968d6SJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
626*3a2968d6SJeremy L Thompson           break;
627*3a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
628*3a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
629*3a2968d6SJeremy L Thompson           code << "      GradAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix
630*3a2968d6SJeremy L Thompson                << ", r_x, r_s" << var_suffix << ");\n";
631*3a2968d6SJeremy L Thompson           break;
632*3a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
633*3a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
634*3a2968d6SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = 1.0;\n";
635*3a2968d6SJeremy L Thompson           break;
636*3a2968d6SJeremy L Thompson           // LCOV_EXCL_START
637*3a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
638*3a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
639*3a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
640*3a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
641*3a2968d6SJeremy L Thompson       }
642*3a2968d6SJeremy L Thompson     }
643*3a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
644*3a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
645*3a2968d6SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
646*3a2968d6SJeremy L Thompson 
647*3a2968d6SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
648*3a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
649*3a2968d6SJeremy L Thompson       // Basis action
650*3a2968d6SJeremy L Thompson       switch (eval_mode) {
651*3a2968d6SJeremy L Thompson         case CEED_EVAL_NONE:
652*3a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
653*3a2968d6SJeremy L Thompson           break;
654*3a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
655*3a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
656*3a2968d6SJeremy L Thompson           break;
657*3a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
658*3a2968d6SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
659*3a2968d6SJeremy L Thompson           break;
660*3a2968d6SJeremy L Thompson           // LCOV_EXCL_START
661*3a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
662*3a2968d6SJeremy L Thompson           break;  // Should not occur
663*3a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
664*3a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
665*3a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
666*3a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
667*3a2968d6SJeremy L Thompson       }
668*3a2968d6SJeremy L Thompson     }
669*3a2968d6SJeremy L Thompson 
670*3a2968d6SJeremy L Thompson   } else if (use_3d_slices) {
6714b3e95d5SJeremy L Thompson     // We treat quadrature points per slice in 3d to save registers
6724b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
6734b3e95d5SJeremy L Thompson     code << "    #pragma unroll\n";
6744b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
6754b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
6764b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
6774b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
6784b3e95d5SJeremy L Thompson 
6794b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
6804b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
6814b3e95d5SJeremy L Thompson       // Basis action
6824b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
6834b3e95d5SJeremy L Thompson       switch (eval_mode) {
6844b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
6854b3e95d5SJeremy L Thompson           bool is_strided;
6864b3e95d5SJeremy L Thompson 
6874b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
6884b3e95d5SJeremy L Thompson 
6894b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
6904b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
6914b3e95d5SJeremy L Thompson           if (is_strided) {
6924b3e95d5SJeremy L Thompson             bool    has_backend_strides;
6934b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
6944b3e95d5SJeremy L Thompson 
6954b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
6964b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
6974b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
6984b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
6994b3e95d5SJeremy L Thompson 
7004b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
7014b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
7024b3e95d5SJeremy L Thompson             }
7034b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
704f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << "," << strides[0] << "," << strides[1] << ","
7054b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
7064b3e95d5SJeremy L Thompson           } else {
7074b3e95d5SJeremy L Thompson             CeedSize                 l_size = 0;
7084b3e95d5SJeremy L Thompson             CeedInt                  comp_stride;
7094b3e95d5SJeremy L Thompson             CeedElemRestriction_Hip *rstr_data;
7104b3e95d5SJeremy L Thompson 
7114b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
7124b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
7134b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
7144b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
7154b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
7164b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
717f815fac9SJeremy L Thompson             code << "      ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
7184b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
7194b3e95d5SJeremy L Thompson           }
7204b3e95d5SJeremy L Thompson           break;
7214b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
7224b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7234b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
7244b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
7254b3e95d5SJeremy L Thompson           code << "      }\n";
7264b3e95d5SJeremy L Thompson           break;
7274b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
7284b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
729f815fac9SJeremy L Thompson           code << "      GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix
730f815fac9SJeremy L Thompson                << ", r_s" << var_suffix << ");\n";
7314b3e95d5SJeremy L Thompson           break;
7324b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7334b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
7344b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
735*3a2968d6SJeremy L Thompson           break;
7364b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
7374b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
7384b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
7394b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
7404b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
7414b3e95d5SJeremy L Thompson       }
7424b3e95d5SJeremy L Thompson     }
7434b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
7444b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
7454b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
7464b3e95d5SJeremy L Thompson 
7474b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
7484b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
7494b3e95d5SJeremy L Thompson       // Basis action
7504b3e95d5SJeremy L Thompson       switch (eval_mode) {
7514b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
7524b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
753*3a2968d6SJeremy L Thompson           break;
7544b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
7554b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
7564b3e95d5SJeremy L Thompson           break;
7574b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
7584b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
7594b3e95d5SJeremy L Thompson           break;
7604b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
7614b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
7624b3e95d5SJeremy L Thompson           break;  // Should not occur
7634b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
7644b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
7654b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
7664b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
7674b3e95d5SJeremy L Thompson       }
7684b3e95d5SJeremy L Thompson     }
7694b3e95d5SJeremy L Thompson   } else {
7704b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
7714b3e95d5SJeremy L Thompson     code << "    {\n";
7724b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
7734b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
7744b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
7754b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
7764b3e95d5SJeremy L Thompson     }
7774b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
7784b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
7794b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
7804b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
7814b3e95d5SJeremy L Thompson     }
7824b3e95d5SJeremy L Thompson   }
7834b3e95d5SJeremy L Thompson 
7844b3e95d5SJeremy L Thompson   // Input and output buffers
7854b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
7864b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
7874b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
7884b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
7894b3e95d5SJeremy L Thompson     code << "      // ------ Input field " << i << "\n";
7904b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
7914b3e95d5SJeremy L Thompson   }
7924b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
7934b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
7944b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
7954b3e95d5SJeremy L Thompson     code << "      // ------ Output field " << i << "\n";
7964b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
7974b3e95d5SJeremy L Thompson   }
7984b3e95d5SJeremy L Thompson 
7994b3e95d5SJeremy L Thompson   // Apply QFunction
8004b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
8014b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
802*3a2968d6SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices) {
8034b3e95d5SJeremy L Thompson     code << "1";
8044b3e95d5SJeremy L Thompson   } else {
8054b3e95d5SJeremy L Thompson     code << "Q_1d";
8064b3e95d5SJeremy L Thompson   }
8074b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
8084b3e95d5SJeremy L Thompson 
809*3a2968d6SJeremy L Thompson   if (is_at_points) {
810*3a2968d6SJeremy L Thompson     // Map back to coefficients
811*3a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
812*3a2968d6SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
813*3a2968d6SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
814*3a2968d6SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
815*3a2968d6SJeremy L Thompson 
816*3a2968d6SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
817*3a2968d6SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
818*3a2968d6SJeremy L Thompson       // Basis action
819*3a2968d6SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
820*3a2968d6SJeremy L Thompson       switch (eval_mode) {
821*3a2968d6SJeremy L Thompson         case CEED_EVAL_NONE: {
822*3a2968d6SJeremy L Thompson           CeedInt             comp_stride;
823*3a2968d6SJeremy L Thompson           CeedElemRestriction elem_rstr;
824*3a2968d6SJeremy L Thompson 
825*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
826*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
827*3a2968d6SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
828*3a2968d6SJeremy L Thompson           code << "      const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n";
829*3a2968d6SJeremy L Thompson           code << "      WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix
830*3a2968d6SJeremy L Thompson                << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]"
831*3a2968d6SJeremy L Thompson                << ", r_s" << var_suffix << ", d" << var_suffix << ");\n";
832*3a2968d6SJeremy L Thompson           break;
833*3a2968d6SJeremy L Thompson         }
834*3a2968d6SJeremy L Thompson         case CEED_EVAL_INTERP:
835*3a2968d6SJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
836*3a2968d6SJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n";
837*3a2968d6SJeremy L Thompson           code << "      }\n";
838*3a2968d6SJeremy L Thompson           code << "      InterpTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
839*3a2968d6SJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
840*3a2968d6SJeremy L Thompson           break;
841*3a2968d6SJeremy L Thompson         case CEED_EVAL_GRAD:
842*3a2968d6SJeremy L Thompson           code << "      if (i >= points.num_per_elem[elem]) {\n";
843*3a2968d6SJeremy L Thompson           code << "        for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim; j++) r_s" << var_suffix << "[j] = 0.0;\n";
844*3a2968d6SJeremy L Thompson           code << "      }\n";
845*3a2968d6SJeremy L Thompson           code << "      GradTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s"
846*3a2968d6SJeremy L Thompson                << var_suffix << ", r_x, r_c" << var_suffix << ");\n";
847*3a2968d6SJeremy L Thompson           break;
848*3a2968d6SJeremy L Thompson           // LCOV_EXCL_START
849*3a2968d6SJeremy L Thompson         case CEED_EVAL_WEIGHT:
850*3a2968d6SJeremy L Thompson           break;  // Should not occur
851*3a2968d6SJeremy L Thompson         case CEED_EVAL_DIV:
852*3a2968d6SJeremy L Thompson         case CEED_EVAL_CURL:
853*3a2968d6SJeremy L Thompson           break;  // TODO: Not implemented
854*3a2968d6SJeremy L Thompson                   // LCOV_EXCL_STOP
855*3a2968d6SJeremy L Thompson       }
856*3a2968d6SJeremy L Thompson     }
857*3a2968d6SJeremy L Thompson   } else if (use_3d_slices) {
8584b3e95d5SJeremy L Thompson     // Copy or apply transpose grad, if needed
859*3a2968d6SJeremy L Thompson     code << "\n      // -- Output fields\n";
8604b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
8614b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
8624b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
8634b3e95d5SJeremy L Thompson 
8644b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
8654b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
8664b3e95d5SJeremy L Thompson       // Basis action
8674b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
8684b3e95d5SJeremy L Thompson       switch (eval_mode) {
8694b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
8704b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
8714b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
8724b3e95d5SJeremy L Thompson           code << "      }\n";
873*3a2968d6SJeremy L Thompson           break;
8744b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
8754b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
8764b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
8774b3e95d5SJeremy L Thompson           code << "      }\n";
8784b3e95d5SJeremy L Thompson           break;
8794b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
880f815fac9SJeremy L Thompson           code << "      GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G"
881f815fac9SJeremy L Thompson                << var_suffix << ", r_q" << var_suffix << ");\n";
8824b3e95d5SJeremy L Thompson           break;
8834b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
8844b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
8854b3e95d5SJeremy L Thompson           break;  // Should not occur
8864b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
8874b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
8884b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
8894b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
8904b3e95d5SJeremy L Thompson       }
8914b3e95d5SJeremy L Thompson     }
8924b3e95d5SJeremy L Thompson   }
8934b3e95d5SJeremy L Thompson   code << "    }\n";
8944b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8954b3e95d5SJeremy L Thompson }
8964b3e95d5SJeremy L Thompson 
8974b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
8989e201c85SYohann // Build single operator kernel
8997d8d0e25Snbeams //------------------------------------------------------------------------------
900eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op) {
901*3a2968d6SJeremy L Thompson   bool                   is_tensor = true, is_at_points = false, use_3d_slices = false;
9027d8d0e25Snbeams   Ceed                   ceed;
903*3a2968d6SJeremy L Thompson   CeedInt                Q_1d, num_input_fields, num_output_fields, dim = 1, max_num_points = 0, coords_comp_stride = 0;
904b7453713SJeremy L Thompson   CeedQFunctionField    *qf_input_fields, *qf_output_fields;
905b7453713SJeremy L Thompson   CeedQFunction_Hip_gen *qf_data;
906b7453713SJeremy L Thompson   CeedQFunction          qf;
907b7453713SJeremy L Thompson   CeedOperatorField     *op_input_fields, *op_output_fields;
908b7453713SJeremy L Thompson   CeedOperator_Hip_gen  *data;
9094b3e95d5SJeremy L Thompson   std::ostringstream     code;
9104b3e95d5SJeremy L Thompson 
9114b3e95d5SJeremy L Thompson   {
9124b3e95d5SJeremy L Thompson     bool is_setup_done;
913b7453713SJeremy L Thompson 
914b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
915b7453713SJeremy L Thompson     if (is_setup_done) return CEED_ERROR_SUCCESS;
9164b3e95d5SJeremy L Thompson   }
917b7453713SJeremy L Thompson 
918b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
919b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
920b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
921b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
922b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
923b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
9247d8d0e25Snbeams 
9254b3e95d5SJeremy L Thompson   // Get operator data
926*3a2968d6SJeremy L Thompson   CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points));
9274b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelData_Hip_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
9284b3e95d5SJeremy L Thompson                                                       qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
9294b3e95d5SJeremy L Thompson   if (dim == 0) dim = 1;
9304b3e95d5SJeremy L Thompson   data->dim = dim;
931*3a2968d6SJeremy L Thompson   if (is_at_points) {
932*3a2968d6SJeremy L Thompson     CeedElemRestriction_Hip *rstr_data;
933*3a2968d6SJeremy L Thompson     CeedElemRestriction      rstr_points = NULL;
9344b3e95d5SJeremy L Thompson 
935*3a2968d6SJeremy L Thompson     CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL));
936*3a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points));
937*3a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride));
938*3a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data));
939*3a2968d6SJeremy L Thompson     data->points.indices = (CeedInt *)rstr_data->d_offsets;
940*3a2968d6SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points));
941*3a2968d6SJeremy L Thompson   }
942*3a2968d6SJeremy L Thompson   if (is_at_points) use_3d_slices = false;
943*3a2968d6SJeremy L Thompson   if (Q_1d == 0) {
944*3a2968d6SJeremy L Thompson     if (is_at_points) Q_1d = max_num_points;
945*3a2968d6SJeremy L Thompson     else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d));
9464b3e95d5SJeremy L Thompson   }
9474b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
9484b3e95d5SJeremy L Thompson 
9490b454692Sjeremylt   // Check for restriction only identity operator
9504b3e95d5SJeremy L Thompson   {
9514b3e95d5SJeremy L Thompson     bool is_identity_qf;
9524b3e95d5SJeremy L Thompson 
9532b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
9540b454692Sjeremylt     if (is_identity_qf) {
9559e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
956b7453713SJeremy L Thompson 
9572b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
9582b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
9596574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
9606574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
9610b454692Sjeremylt     }
9624b3e95d5SJeremy L Thompson   }
963b2165e7aSSebastian Grimberg 
964b2165e7aSSebastian Grimberg   // Load basis source files
9654b3e95d5SJeremy L Thompson   // TODO: Add non-tensor, AtPoints
9669c25dd66SJeremy L Thompson   code << "// Tensor basis source\n";
9679c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-templates.h>\n\n";
968*3a2968d6SJeremy L Thompson   code << "// AtPoints basis source\n";
969*3a2968d6SJeremy L Thompson   code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-at-points-templates.h>\n\n";
9709c25dd66SJeremy L Thompson   code << "// CodeGen operator source\n";
9719c25dd66SJeremy L Thompson   code << "#include <ceed/jit-source/hip/hip-gen-templates.h>\n\n";
9727d8d0e25Snbeams 
9734b3e95d5SJeremy L Thompson   // Get QFunction name
9744b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
9754b3e95d5SJeremy L Thompson   std::string operator_name;
9764b3e95d5SJeremy L Thompson 
97709095acaSJeremy L Thompson   operator_name = "CeedKernelHipGenOperator_" + qfunction_name;
9787d8d0e25Snbeams 
9799e201c85SYohann   // Define CEED_Q_VLA
9809e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
981*3a2968d6SJeremy L Thompson   if (dim != 3 || is_at_points || use_3d_slices) {
9829e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
9839e201c85SYohann   } else {
9849e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
9859e201c85SYohann   }
9869e201c85SYohann 
9874b3e95d5SJeremy L Thompson   // Add user QFunction source
9884b3e95d5SJeremy L Thompson   {
9899c25dd66SJeremy L Thompson     const char *source_path;
9904b3e95d5SJeremy L Thompson 
9919c25dd66SJeremy L Thompson     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
9929c25dd66SJeremy L Thompson     CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/hip/gen backend requires QFunction source code file");
9939c25dd66SJeremy L Thompson 
9949c25dd66SJeremy L Thompson     code << "// User QFunction source\n";
9959c25dd66SJeremy L Thompson     code << "#include \"" << source_path << "\"\n\n";
9964b3e95d5SJeremy L Thompson   }
9977d8d0e25Snbeams 
9987d8d0e25Snbeams   // Setup
9997d8d0e25Snbeams   code << "\n// -----------------------------------------------------------------------------\n";
10004b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
10014b3e95d5SJeremy L Thompson   code << "// \n";
10024b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
10034b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
10044b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
1005*3a2968d6SJeremy L Thompson   code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficents register\n";
10064b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice vector register\n";
10074b3e95d5SJeremy L Thompson   code << "// \n";
10084b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
10094b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
10104b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
1011b3e1519bSnbeams   code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
10122b730f8bSJeremy L Thompson   code << "__global__ void " << operator_name
1013*3a2968d6SJeremy 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";
10144b3e95d5SJeremy L Thompson 
10154b3e95d5SJeremy L Thompson   // Scratch buffers
10169e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
10174b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
10184b3e95d5SJeremy L Thompson 
10192b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
10209e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
10214b3e95d5SJeremy L Thompson       code << "  const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n";
10227d8d0e25Snbeams     }
10237d8d0e25Snbeams   }
10249e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
10254b3e95d5SJeremy L Thompson     code << "  CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n";
10267d8d0e25Snbeams   }
10277d8d0e25Snbeams 
10289e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
10299e201c85SYohann   code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
1030*3a2968d6SJeremy L Thompson   if (is_at_points) {
1031*3a2968d6SJeremy L Thompson     code << "  const CeedInt max_num_points = " << max_num_points << ";\n";
1032*3a2968d6SJeremy L Thompson     code << "  const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n";
1033*3a2968d6SJeremy L Thompson   }
10347d8d0e25Snbeams 
10354b3e95d5SJeremy L Thompson   // Shared data
10364b3e95d5SJeremy L Thompson   code << "  extern __shared__ CeedScalar slice[];\n";
10379e201c85SYohann   code << "  SharedData_Hip data;\n";
10389e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
10399e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
10409e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
10419e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
10429e201c85SYohann   code << "  data.slice = slice + data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n";
10437d8d0e25Snbeams 
10447d8d0e25Snbeams   // Initialize constants, and matrices B and G
10454b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
10469e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1047*3a2968d6SJeremy L Thompson     CeedCallBackend(
1048*3a2968d6SJeremy L Thompson         CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_input_fields[i], qf_input_fields[i], Q_1d, true, is_at_points, use_3d_slices));
10497d8d0e25Snbeams   }
10504b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
10519e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1052*3a2968d6SJeremy L Thompson     CeedCallBackend(
1053*3a2968d6SJeremy L Thompson         CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_at_points, use_3d_slices));
10544b3e95d5SJeremy L Thompson   }
10557d8d0e25Snbeams 
10564b3e95d5SJeremy L Thompson   // Loop over all elements
10574b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
10587d8d0e25Snbeams   code << "  __syncthreads();\n";
10599e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
10604b3e95d5SJeremy L Thompson 
1061e93651e5SJeremy L Thompson   // -- Compute minimum buffer space needed
1062*3a2968d6SJeremy L Thompson   CeedInt max_rstr_buffer_size = 1;
1063e93651e5SJeremy L Thompson 
1064e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1065e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1066e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1067e93651e5SJeremy L Thompson 
1068e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
1069e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1070e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1071e93651e5SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size);
1072681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1073e93651e5SJeremy L Thompson   }
1074e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
1075e93651e5SJeremy L Thompson     CeedInt             num_comp, elem_size;
1076e93651e5SJeremy L Thompson     CeedElemRestriction elem_rstr;
1077e93651e5SJeremy L Thompson 
1078e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr));
1079e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
1080e93651e5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
1081e93651e5SJeremy L Thompson     max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size);
1082681d0ea7SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr));
1083e93651e5SJeremy L Thompson   }
1084e93651e5SJeremy L Thompson   code << "    // Scratch restriction buffer space\n";
1085e93651e5SJeremy L Thompson   code << "    CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n";
1086e93651e5SJeremy L Thompson 
1087e93651e5SJeremy L Thompson   // -- Determine best input field processing order
1088e93651e5SJeremy L Thompson   CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX];
1089e93651e5SJeremy L Thompson 
1090e93651e5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
1091e93651e5SJeremy L Thompson     field_rstr_in_buffer[i] = -1;
1092e93651e5SJeremy L Thompson     input_field_order[i]    = -1;
1093e93651e5SJeremy L Thompson   }
1094e93651e5SJeremy L Thompson   {
1095e93651e5SJeremy L Thompson     bool    is_ordered[CEED_FIELD_MAX];
1096e93651e5SJeremy L Thompson     CeedInt curr_index = 0;
1097e93651e5SJeremy L Thompson 
1098e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false;
1099e93651e5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
1100e93651e5SJeremy L Thompson       CeedVector          vec_i;
1101e93651e5SJeremy L Thompson       CeedElemRestriction rstr_i;
1102e93651e5SJeremy L Thompson 
1103e93651e5SJeremy L Thompson       if (is_ordered[i]) continue;
1104e93651e5SJeremy L Thompson       field_rstr_in_buffer[i]       = i;
1105e93651e5SJeremy L Thompson       is_ordered[i]                 = true;
1106e93651e5SJeremy L Thompson       input_field_order[curr_index] = i;
1107e93651e5SJeremy L Thompson       curr_index++;
1108034f99fdSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i));
1109e93651e5SJeremy L Thompson       if (vec_i == CEED_VECTOR_NONE) continue;  // CEED_EVAL_WEIGHT
1110e93651e5SJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i));
1111e93651e5SJeremy L Thompson       for (CeedInt j = i + 1; j < num_input_fields; j++) {
1112e93651e5SJeremy L Thompson         CeedVector          vec_j;
1113e93651e5SJeremy L Thompson         CeedElemRestriction rstr_j;
1114e93651e5SJeremy L Thompson 
1115e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j));
1116e93651e5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j));
1117e93651e5SJeremy L Thompson         if (rstr_i == rstr_j && vec_i == vec_j) {
1118e93651e5SJeremy L Thompson           field_rstr_in_buffer[j]       = i;
1119e93651e5SJeremy L Thompson           is_ordered[j]                 = true;
1120e93651e5SJeremy L Thompson           input_field_order[curr_index] = j;
1121e93651e5SJeremy L Thompson           curr_index++;
1122e93651e5SJeremy L Thompson         }
1123*3a2968d6SJeremy L Thompson         CeedCallBackend(CeedVectorDestroy(&vec_j));
1124*3a2968d6SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j));
1125e93651e5SJeremy L Thompson       }
1126*3a2968d6SJeremy L Thompson       CeedCallBackend(CeedVectorDestroy(&vec_i));
1127*3a2968d6SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i));
1128e93651e5SJeremy L Thompson     }
1129e93651e5SJeremy L Thompson   }
1130e93651e5SJeremy L Thompson 
11314b3e95d5SJeremy L Thompson   // -- Input restriction and basis
1132*3a2968d6SJeremy L Thompson   code << "\n    // -- Input field restrictions and basis actions\n";
11339e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1134e93651e5SJeremy L Thompson     CeedInt f = input_field_order[i];
1135e93651e5SJeremy L Thompson 
1136e93651e5SJeremy L Thompson     code << "    // ---- Input field " << f << "\n";
11377d8d0e25Snbeams 
11384b3e95d5SJeremy L Thompson     // ---- Restriction
1139e93651e5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], Q_1d,
1140*3a2968d6SJeremy L Thompson                                                                true, is_at_points, use_3d_slices));
1141b7453713SJeremy L Thompson 
11424b3e95d5SJeremy L Thompson     // ---- Basis action
1143*3a2968d6SJeremy L Thompson     CeedCallBackend(
1144*3a2968d6SJeremy L Thompson         CeedOperatorBuildKernelBasis_Hip_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, is_at_points, use_3d_slices));
11457d8d0e25Snbeams   }
11467d8d0e25Snbeams 
11474b3e95d5SJeremy L Thompson   // -- Q function
1148*3a2968d6SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Hip_gen(code, data, dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields,
1149*3a2968d6SJeremy L Thompson                                                            num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d, is_at_points,
1150*3a2968d6SJeremy L Thompson                                                            use_3d_slices));
11517d8d0e25Snbeams 
11524b3e95d5SJeremy L Thompson   // -- Output basis and restriction
11534b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
11549e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
11554b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
1156b7453713SJeremy L Thompson 
11574b3e95d5SJeremy L Thompson     // ---- Basis action
1158*3a2968d6SJeremy L Thompson     CeedCallBackend(
1159*3a2968d6SJeremy L Thompson         CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_at_points, use_3d_slices));
11607d8d0e25Snbeams 
11614b3e95d5SJeremy L Thompson     // ---- Restriction
1162*3a2968d6SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false,
1163*3a2968d6SJeremy L Thompson                                                                is_at_points, use_3d_slices));
11647d8d0e25Snbeams   }
11657d8d0e25Snbeams 
11664b3e95d5SJeremy L Thompson   // Close loop and function
11677d8d0e25Snbeams   code << "  }\n";
11687d8d0e25Snbeams   code << "}\n";
11697d8d0e25Snbeams   code << "// -----------------------------------------------------------------------------\n\n";
11707d8d0e25Snbeams 
1171539ec17dSJeremy L Thompson   CeedInt block_sizes[3] = {0, 0, 0};
11729e201c85SYohann   CeedInt num_elem;
1173b7453713SJeremy L Thompson 
1174*3a2968d6SJeremy L Thompson   // Compile
11752b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
11762b730f8bSJeremy L Thompson   CeedCallBackend(BlockGridCalculate_Hip_gen(dim, num_elem, data->max_P_1d, Q_1d, block_sizes));
1177eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE",
11782b730f8bSJeremy L Thompson                                   block_sizes[0] * block_sizes[1] * block_sizes[2]));
1179eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op));
11802b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
11819bc66399SJeremy L Thompson   CeedCallBackend(CeedDestroy(&ceed));
1182c11e12f4SJeremy L Thompson   CeedCallBackend(CeedQFunctionDestroy(&qf));
1183e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
11847d8d0e25Snbeams }
11852a86cc9dSSebastian Grimberg 
11867d8d0e25Snbeams //------------------------------------------------------------------------------
1187