xref: /libCEED/rust/libceed-sys/c-src/backends/hip-gen/ceed-hip-gen-operator-build.cpp (revision 4b3e95d5da7d4d1a3e8f0acb9d6f6cc36918381b)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
37d8d0e25Snbeams //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
57d8d0e25Snbeams //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
73d576824SJeremy L Thompson 
87d8d0e25Snbeams #define CEED_DEBUG_COLOR 12
97d8d0e25Snbeams 
1049aac155SJeremy L Thompson #include <ceed.h>
11ec3da8bcSJed Brown #include <ceed/backend.h>
129e201c85SYohann #include <ceed/jit-tools.h>
132b730f8bSJeremy L Thompson 
147d8d0e25Snbeams #include <iostream>
157d8d0e25Snbeams #include <sstream>
162b730f8bSJeremy L Thompson #include <string>
172b730f8bSJeremy L Thompson 
180d0321e0SJeremy L Thompson #include "../hip-ref/ceed-hip-ref.h"
197d8d0e25Snbeams #include "../hip-shared/ceed-hip-shared.h"
20b2165e7aSSebastian Grimberg #include "../hip/ceed-hip-common.h"
217d8d0e25Snbeams #include "../hip/ceed-hip-compile.h"
222b730f8bSJeremy L Thompson #include "ceed-hip-gen.h"
23b3e1519bSnbeams 
24b3e1519bSnbeams //------------------------------------------------------------------------------
25b3e1519bSnbeams // Calculate the block size used for launching the operator kernel
26b3e1519bSnbeams //------------------------------------------------------------------------------
272b730f8bSJeremy L Thompson extern "C" int BlockGridCalculate_Hip_gen(const CeedInt dim, const CeedInt num_elem, const CeedInt P_1d, const CeedInt Q_1d, CeedInt *block_sizes) {
289e201c85SYohann   const CeedInt thread1d = CeedIntMax(Q_1d, P_1d);
29b3e1519bSnbeams   if (dim == 1) {
309e201c85SYohann     CeedInt elems_per_block = 64 * thread1d > 256 ? 256 / thread1d : 64;
31b7453713SJeremy L Thompson 
329e201c85SYohann     elems_per_block = elems_per_block > 0 ? elems_per_block : 1;
33b3e1519bSnbeams     block_sizes[0]  = thread1d;
34b3e1519bSnbeams     block_sizes[1]  = 1;
359e201c85SYohann     block_sizes[2]  = elems_per_block;
36b3e1519bSnbeams   } else if (dim == 2) {
379e201c85SYohann     const CeedInt elems_per_block = thread1d < 4 ? 16 : 2;
38b7453713SJeremy L Thompson 
39b3e1519bSnbeams     block_sizes[0] = thread1d;
40b3e1519bSnbeams     block_sizes[1] = thread1d;
419e201c85SYohann     block_sizes[2] = elems_per_block;
42b3e1519bSnbeams   } else if (dim == 3) {
439e201c85SYohann     const CeedInt elems_per_block = thread1d < 6 ? 4 : (thread1d < 8 ? 2 : 1);
44b7453713SJeremy L Thompson 
45b3e1519bSnbeams     block_sizes[0] = thread1d;
46b3e1519bSnbeams     block_sizes[1] = thread1d;
479e201c85SYohann     block_sizes[2] = elems_per_block;
48b3e1519bSnbeams   }
49b3e1519bSnbeams   return CEED_ERROR_SUCCESS;
50b3e1519bSnbeams }
51b3e1519bSnbeams 
527d8d0e25Snbeams //------------------------------------------------------------------------------
53*4b3e95d5SJeremy L Thompson // Determine type of operator
54*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
55*4b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Hip_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields,
56*4b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields,
57*4b3e95d5SJeremy L Thompson                                                CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor,
58*4b3e95d5SJeremy L Thompson                                                bool *use_3d_slices) {
59*4b3e95d5SJeremy L Thompson   // Find dim, P_1d, Q_1d
60*4b3e95d5SJeremy L Thompson   *max_P_1d  = 0;
61*4b3e95d5SJeremy L Thompson   *Q_1d      = 0;
62*4b3e95d5SJeremy L Thompson   *dim       = 0;
63*4b3e95d5SJeremy L Thompson   *is_tensor = true;
64*4b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
65*4b3e95d5SJeremy L Thompson     CeedBasis basis;
66*4b3e95d5SJeremy L Thompson 
67*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
68*4b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
69*4b3e95d5SJeremy L Thompson       bool    is_field_tensor;
70*4b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
71*4b3e95d5SJeremy L Thompson 
72*4b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
73*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
74*4b3e95d5SJeremy L Thompson       CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
75*4b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
76*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
77*4b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
78*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
79*4b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
80*4b3e95d5SJeremy L Thompson       *dim = field_dim;
81*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
82*4b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
83*4b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
84*4b3e95d5SJeremy L Thompson     }
85*4b3e95d5SJeremy L Thompson   }
86*4b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
87*4b3e95d5SJeremy L Thompson     CeedBasis basis;
88*4b3e95d5SJeremy L Thompson 
89*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
90*4b3e95d5SJeremy L Thompson     if (basis != CEED_BASIS_NONE) {
91*4b3e95d5SJeremy L Thompson       bool    is_field_tensor;
92*4b3e95d5SJeremy L Thompson       CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0;
93*4b3e95d5SJeremy L Thompson 
94*4b3e95d5SJeremy L Thompson       // Collect dim, P_1d, and Q_1d
95*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor));
96*4b3e95d5SJeremy L Thompson       CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis");
97*4b3e95d5SJeremy L Thompson       *is_tensor = *is_tensor && is_field_tensor;
98*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d));
99*4b3e95d5SJeremy L Thompson       *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d);
100*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetDimension(basis, &field_dim));
101*4b3e95d5SJeremy L Thompson       CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
102*4b3e95d5SJeremy L Thompson       *dim = field_dim;
103*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d));
104*4b3e95d5SJeremy L Thompson       CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible");
105*4b3e95d5SJeremy L Thompson       *Q_1d = field_Q_1d;
106*4b3e95d5SJeremy L Thompson     }
107*4b3e95d5SJeremy L Thompson   }
108*4b3e95d5SJeremy L Thompson 
109*4b3e95d5SJeremy L Thompson   // Only use 3D collocated gradient parallelization strategy when gradient is computed
110*4b3e95d5SJeremy L Thompson   *use_3d_slices = false;
111*4b3e95d5SJeremy L Thompson   if (*dim == 3) {
112*4b3e95d5SJeremy L Thompson     bool was_grad_found = false;
113*4b3e95d5SJeremy L Thompson 
114*4b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
115*4b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
116*4b3e95d5SJeremy L Thompson 
117*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
118*4b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
119*4b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
120*4b3e95d5SJeremy L Thompson         CeedBasis             basis;
121*4b3e95d5SJeremy L Thompson 
122*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis));
123*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
124*4b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
125*4b3e95d5SJeremy L Thompson         was_grad_found = true;
126*4b3e95d5SJeremy L Thompson       }
127*4b3e95d5SJeremy L Thompson     }
128*4b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
129*4b3e95d5SJeremy L Thompson       CeedEvalMode eval_mode;
130*4b3e95d5SJeremy L Thompson 
131*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
132*4b3e95d5SJeremy L Thompson       if (eval_mode == CEED_EVAL_GRAD) {
133*4b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
134*4b3e95d5SJeremy L Thompson         CeedBasis             basis;
135*4b3e95d5SJeremy L Thompson 
136*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis));
137*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
138*4b3e95d5SJeremy L Thompson         *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true);
139*4b3e95d5SJeremy L Thompson         was_grad_found = true;
140*4b3e95d5SJeremy L Thompson       }
141*4b3e95d5SJeremy L Thompson     }
142*4b3e95d5SJeremy L Thompson   }
143*4b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
144*4b3e95d5SJeremy L Thompson }
145*4b3e95d5SJeremy L Thompson 
146*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
147*4b3e95d5SJeremy L Thompson // Setup fields
148*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
149*4b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedOperatorField op_field,
150*4b3e95d5SJeremy L Thompson                                                     CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool use_3d_slices) {
151*4b3e95d5SJeremy L Thompson   std::string           var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
152*4b3e95d5SJeremy L Thompson   std::string           P_name = "P_1d" + var_suffix, Q_name = "Q_1d";
153*4b3e95d5SJeremy L Thompson   std::string           option_name = (is_input ? "inputs" : "outputs");
154*4b3e95d5SJeremy L Thompson   CeedEvalMode          eval_mode   = CEED_EVAL_NONE;
155*4b3e95d5SJeremy L Thompson   CeedInt               elem_size = 0, num_comp = 0, P_1d = 0;
156*4b3e95d5SJeremy L Thompson   CeedElemRestriction   elem_rstr;
157*4b3e95d5SJeremy L Thompson   CeedBasis_Hip_shared *basis_data;
158*4b3e95d5SJeremy L Thompson   CeedBasis             basis;
159*4b3e95d5SJeremy L Thompson 
160*4b3e95d5SJeremy L Thompson   code << "  // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n";
161*4b3e95d5SJeremy L Thompson 
162*4b3e95d5SJeremy L Thompson   // Get field data
163*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
164*4b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
165*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
166*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
167*4b3e95d5SJeremy L Thompson   }
168*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
169*4b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
170*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetData(basis, &basis_data));
171*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
172*4b3e95d5SJeremy L Thompson   }
173*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
174*4b3e95d5SJeremy L Thompson 
175*4b3e95d5SJeremy L Thompson   // Set field constants
176*4b3e95d5SJeremy L Thompson   if (eval_mode != CEED_EVAL_WEIGHT) {
177*4b3e95d5SJeremy L Thompson     code << "  const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n";
178*4b3e95d5SJeremy L Thompson     code << "  const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n";
179*4b3e95d5SJeremy L Thompson   }
180*4b3e95d5SJeremy L Thompson 
181*4b3e95d5SJeremy L Thompson   // Load basis data
182*4b3e95d5SJeremy L Thompson   code << "  // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
183*4b3e95d5SJeremy L Thompson   switch (eval_mode) {
184*4b3e95d5SJeremy L Thompson     case CEED_EVAL_NONE:
185*4b3e95d5SJeremy L Thompson       break;
186*4b3e95d5SJeremy L Thompson     case CEED_EVAL_INTERP:
187*4b3e95d5SJeremy L Thompson       if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
188*4b3e95d5SJeremy L Thompson       else data->B.outputs[i] = basis_data->d_interp_1d;
189*4b3e95d5SJeremy L Thompson       code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n";
190*4b3e95d5SJeremy L Thompson       code << "  loadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
191*4b3e95d5SJeremy L Thompson       break;
192*4b3e95d5SJeremy L Thompson     case CEED_EVAL_GRAD:
193*4b3e95d5SJeremy L Thompson       if (is_input) data->B.inputs[i] = basis_data->d_interp_1d;
194*4b3e95d5SJeremy L Thompson       else data->B.outputs[i] = basis_data->d_interp_1d;
195*4b3e95d5SJeremy L Thompson       code << "  __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n";
196*4b3e95d5SJeremy L Thompson       code << "  loadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n";
197*4b3e95d5SJeremy L Thompson       if (use_3d_slices) {
198*4b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d;
199*4b3e95d5SJeremy L Thompson         else data->G.outputs[i] = basis_data->d_collo_grad_1d;
200*4b3e95d5SJeremy L Thompson         code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n";
201*4b3e95d5SJeremy L Thompson         code << "  loadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
202*4b3e95d5SJeremy L Thompson       } else {
203*4b3e95d5SJeremy L Thompson         bool has_collo_grad = basis_data->d_collo_grad_1d;
204*4b3e95d5SJeremy L Thompson 
205*4b3e95d5SJeremy L Thompson         if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
206*4b3e95d5SJeremy L Thompson         else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d;
207*4b3e95d5SJeremy L Thompson         if (has_collo_grad) {
208*4b3e95d5SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n";
209*4b3e95d5SJeremy L Thompson           code << "  loadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
210*4b3e95d5SJeremy L Thompson         } else {
211*4b3e95d5SJeremy L Thompson           code << "  __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * P_1d << "];\n";
212*4b3e95d5SJeremy L Thompson           code << "  loadMatrix<" << P_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n";
213*4b3e95d5SJeremy L Thompson         }
214*4b3e95d5SJeremy L Thompson       }
215*4b3e95d5SJeremy L Thompson       break;
216*4b3e95d5SJeremy L Thompson     case CEED_EVAL_WEIGHT:
217*4b3e95d5SJeremy L Thompson       break;  // No action
218*4b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
219*4b3e95d5SJeremy L Thompson     case CEED_EVAL_DIV:
220*4b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
221*4b3e95d5SJeremy L Thompson     case CEED_EVAL_CURL:
222*4b3e95d5SJeremy L Thompson       break;  // TODO: Not implemented
223*4b3e95d5SJeremy L Thompson               // LCOV_EXCL_STOP
224*4b3e95d5SJeremy L Thompson   }
225*4b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
226*4b3e95d5SJeremy L Thompson }
227*4b3e95d5SJeremy L Thompson 
228*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
229*4b3e95d5SJeremy L Thompson // Restriction
230*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
231*4b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
232*4b3e95d5SJeremy L Thompson                                                       CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input,
233*4b3e95d5SJeremy L Thompson                                                       bool use_3d_slices) {
234*4b3e95d5SJeremy L Thompson   std::string              var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
235*4b3e95d5SJeremy L Thompson   std::string              P_name     = "P_1d" + var_suffix;
236*4b3e95d5SJeremy L Thompson   CeedEvalMode             eval_mode  = CEED_EVAL_NONE;
237*4b3e95d5SJeremy L Thompson   CeedInt                  elem_size = 0, num_comp = 0, P_1d = 0;
238*4b3e95d5SJeremy L Thompson   CeedSize                 l_size;
239*4b3e95d5SJeremy L Thompson   CeedElemRestriction_Hip *rstr_data;
240*4b3e95d5SJeremy L Thompson   CeedElemRestriction      elem_rstr;
241*4b3e95d5SJeremy L Thompson   CeedBasis                basis;
242*4b3e95d5SJeremy L Thompson 
243*4b3e95d5SJeremy L Thompson   // Get field data
244*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
245*4b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
246*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
247*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
248*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
249*4b3e95d5SJeremy L Thompson   }
250*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
251*4b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
252*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
253*4b3e95d5SJeremy L Thompson   }
254*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
255*4b3e95d5SJeremy L Thompson 
256*4b3e95d5SJeremy L Thompson   // Restriction
257*4b3e95d5SJeremy L Thompson   if (is_input) {
258*4b3e95d5SJeremy L Thompson     // Input
259*4b3e95d5SJeremy L Thompson     if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices)) {
260*4b3e95d5SJeremy L Thompson       bool is_strided;
261*4b3e95d5SJeremy L Thompson 
262*4b3e95d5SJeremy L Thompson       code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
263*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
264*4b3e95d5SJeremy L Thompson       if (!is_strided) {
265*4b3e95d5SJeremy L Thompson         CeedInt comp_stride;
266*4b3e95d5SJeremy L Thompson 
267*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
268*4b3e95d5SJeremy L Thompson         code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
269*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
270*4b3e95d5SJeremy L Thompson         code << "    // CompStride: " << comp_stride << "\n";
271*4b3e95d5SJeremy L Thompson         data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
272*4b3e95d5SJeremy L Thompson         code << "    readDofsOffset" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" << var_suffix
273*4b3e95d5SJeremy L Thompson              << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n";
274*4b3e95d5SJeremy L Thompson       } else {
275*4b3e95d5SJeremy L Thompson         bool    has_backend_strides;
276*4b3e95d5SJeremy L Thompson         CeedInt num_elem;
277*4b3e95d5SJeremy L Thompson 
278*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
279*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
280*4b3e95d5SJeremy L Thompson         CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
281*4b3e95d5SJeremy L Thompson 
282*4b3e95d5SJeremy L Thompson         if (!has_backend_strides) {
283*4b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
284*4b3e95d5SJeremy L Thompson         }
285*4b3e95d5SJeremy L Thompson         code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
286*4b3e95d5SJeremy L Thompson         code << "    readDofsStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << ","
287*4b3e95d5SJeremy L Thompson              << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n";
288*4b3e95d5SJeremy L Thompson       }
289*4b3e95d5SJeremy L Thompson     }
290*4b3e95d5SJeremy L Thompson   } else {
291*4b3e95d5SJeremy L Thompson     // Output
292*4b3e95d5SJeremy L Thompson     bool is_strided;
293*4b3e95d5SJeremy L Thompson 
294*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
295*4b3e95d5SJeremy L Thompson     if (!is_strided) {
296*4b3e95d5SJeremy L Thompson       CeedInt comp_stride;
297*4b3e95d5SJeremy L Thompson 
298*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
299*4b3e95d5SJeremy L Thompson       code << "    const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
300*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
301*4b3e95d5SJeremy L Thompson       code << "    // CompStride: " << comp_stride << "\n";
302*4b3e95d5SJeremy L Thompson       data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets;
303*4b3e95d5SJeremy L Thompson       code << "    writeDofsOffset" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" << var_suffix
304*4b3e95d5SJeremy L Thompson            << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n";
305*4b3e95d5SJeremy L Thompson     } else {
306*4b3e95d5SJeremy L Thompson       bool    has_backend_strides;
307*4b3e95d5SJeremy L Thompson       CeedInt num_elem;
308*4b3e95d5SJeremy L Thompson 
309*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
310*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
311*4b3e95d5SJeremy L Thompson       CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
312*4b3e95d5SJeremy L Thompson 
313*4b3e95d5SJeremy L Thompson       if (!has_backend_strides) {
314*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
315*4b3e95d5SJeremy L Thompson       }
316*4b3e95d5SJeremy L Thompson       code << "    // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
317*4b3e95d5SJeremy L Thompson       code << "    writeDofsStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << ","
318*4b3e95d5SJeremy L Thompson            << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n";
319*4b3e95d5SJeremy L Thompson     }
320*4b3e95d5SJeremy L Thompson   }
321*4b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
322*4b3e95d5SJeremy L Thompson }
323*4b3e95d5SJeremy L Thompson 
324*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
325*4b3e95d5SJeremy L Thompson // Basis
326*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
327*4b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim,
328*4b3e95d5SJeremy L Thompson                                                 CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input,
329*4b3e95d5SJeremy L Thompson                                                 bool use_3d_slices) {
330*4b3e95d5SJeremy L Thompson   std::string         var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i);
331*4b3e95d5SJeremy L Thompson   std::string         P_name = "P_1d" + var_suffix, Q_name = "Q_1d";
332*4b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
333*4b3e95d5SJeremy L Thompson   CeedInt             elem_size = 0, num_comp = 0, P_1d = 0;
334*4b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
335*4b3e95d5SJeremy L Thompson   CeedBasis           basis;
336*4b3e95d5SJeremy L Thompson 
337*4b3e95d5SJeremy L Thompson   // Get field data
338*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr));
339*4b3e95d5SJeremy L Thompson   if (elem_rstr != CEED_ELEMRESTRICTION_NONE) {
340*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
341*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp));
342*4b3e95d5SJeremy L Thompson   }
343*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis));
344*4b3e95d5SJeremy L Thompson   if (basis != CEED_BASIS_NONE) {
345*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
346*4b3e95d5SJeremy L Thompson   }
347*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode));
348*4b3e95d5SJeremy L Thompson 
349*4b3e95d5SJeremy L Thompson   // Basis
350*4b3e95d5SJeremy L Thompson   code << "    // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
351*4b3e95d5SJeremy L Thompson   if (is_input) {
352*4b3e95d5SJeremy L Thompson     switch (eval_mode) {
353*4b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
354*4b3e95d5SJeremy L Thompson         if (!use_3d_slices) {
355*4b3e95d5SJeremy L Thompson           code << "    CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n";
356*4b3e95d5SJeremy L Thompson         }
357*4b3e95d5SJeremy L Thompson         break;
358*4b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
359*4b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
360*4b3e95d5SJeremy L Thompson         code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
361*4b3e95d5SJeremy L Thompson              << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
362*4b3e95d5SJeremy L Thompson         break;
363*4b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
364*4b3e95d5SJeremy L Thompson         if (use_3d_slices) {
365*4b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
366*4b3e95d5SJeremy L Thompson           code << "    Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name
367*4b3e95d5SJeremy L Thompson                << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n";
368*4b3e95d5SJeremy L Thompson         } else {
369*4b3e95d5SJeremy L Thompson           code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n";
370*4b3e95d5SJeremy L Thompson           code << "    Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" << var_suffix
371*4b3e95d5SJeremy L Thompson                << ", P_1d" << var_suffix << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q"
372*4b3e95d5SJeremy L Thompson                << var_suffix << ");\n";
373*4b3e95d5SJeremy L Thompson         }
374*4b3e95d5SJeremy L Thompson         break;
375*4b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT: {
376*4b3e95d5SJeremy L Thompson         CeedBasis_Hip_shared *basis_data;
377*4b3e95d5SJeremy L Thompson 
378*4b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[" << Q_name << "];\n";
379*4b3e95d5SJeremy L Thompson         CeedCallBackend(CeedBasisGetData(basis, &basis_data));
380*4b3e95d5SJeremy L Thompson         data->W = basis_data->d_q_weight_1d;
381*4b3e95d5SJeremy L Thompson         code << "    Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n";
382*4b3e95d5SJeremy L Thompson         break;
383*4b3e95d5SJeremy L Thompson       }
384*4b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
385*4b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
386*4b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
387*4b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
388*4b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
389*4b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
390*4b3e95d5SJeremy L Thompson     }
391*4b3e95d5SJeremy L Thompson   } else {
392*4b3e95d5SJeremy L Thompson     switch (eval_mode) {
393*4b3e95d5SJeremy L Thompson       case CEED_EVAL_NONE:
394*4b3e95d5SJeremy L Thompson         code << "    CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n";
395*4b3e95d5SJeremy L Thompson         break;  // No action
396*4b3e95d5SJeremy L Thompson       case CEED_EVAL_INTERP:
397*4b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
398*4b3e95d5SJeremy L Thompson         code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
399*4b3e95d5SJeremy L Thompson              << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
400*4b3e95d5SJeremy L Thompson         break;
401*4b3e95d5SJeremy L Thompson       case CEED_EVAL_GRAD:
402*4b3e95d5SJeremy L Thompson         code << "    CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n";
403*4b3e95d5SJeremy L Thompson         if (use_3d_slices) {
404*4b3e95d5SJeremy L Thompson           code << "    InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name
405*4b3e95d5SJeremy L Thompson                << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n";
406*4b3e95d5SJeremy L Thompson         } else {
407*4b3e95d5SJeremy L Thompson           code << "    GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp"
408*4b3e95d5SJeremy L Thompson                << var_suffix << ", " << P_name << "," << Q_name << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix
409*4b3e95d5SJeremy L Thompson                << ", r_e" << var_suffix << ");\n";
410*4b3e95d5SJeremy L Thompson         }
411*4b3e95d5SJeremy L Thompson         break;
412*4b3e95d5SJeremy L Thompson       // LCOV_EXCL_START
413*4b3e95d5SJeremy L Thompson       case CEED_EVAL_WEIGHT:
414*4b3e95d5SJeremy L Thompson         break;  // Should not occur
415*4b3e95d5SJeremy L Thompson       case CEED_EVAL_DIV:
416*4b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
417*4b3e95d5SJeremy L Thompson       case CEED_EVAL_CURL:
418*4b3e95d5SJeremy L Thompson         break;  // TODO: Not implemented
419*4b3e95d5SJeremy L Thompson                 // LCOV_EXCL_STOP
420*4b3e95d5SJeremy L Thompson     }
421*4b3e95d5SJeremy L Thompson   }
422*4b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
423*4b3e95d5SJeremy L Thompson }
424*4b3e95d5SJeremy L Thompson 
425*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
426*4b3e95d5SJeremy L Thompson // QFunction
427*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
428*4b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt dim, CeedInt num_input_fields,
429*4b3e95d5SJeremy L Thompson                                                     CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields,
430*4b3e95d5SJeremy L Thompson                                                     CeedInt num_output_fields, CeedOperatorField *op_output_fields,
431*4b3e95d5SJeremy L Thompson                                                     CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d,
432*4b3e95d5SJeremy L Thompson                                                     bool use_3d_slices) {
433*4b3e95d5SJeremy L Thompson   std::string         Q_name    = "Q_1d";
434*4b3e95d5SJeremy L Thompson   CeedEvalMode        eval_mode = CEED_EVAL_NONE;
435*4b3e95d5SJeremy L Thompson   CeedElemRestriction elem_rstr;
436*4b3e95d5SJeremy L Thompson 
437*4b3e95d5SJeremy L Thompson   // Setup output arays
438*4b3e95d5SJeremy L Thompson   code << "\n    // -- Output field setup\n";
439*4b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
440*4b3e95d5SJeremy L Thompson     std::string var_suffix = "_out_" + std::to_string(i);
441*4b3e95d5SJeremy L Thompson 
442*4b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
443*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
444*4b3e95d5SJeremy L Thompson     if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) {
445*4b3e95d5SJeremy L Thompson       code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
446*4b3e95d5SJeremy L Thompson     }
447*4b3e95d5SJeremy L Thompson     if (eval_mode == CEED_EVAL_GRAD) {
448*4b3e95d5SJeremy L Thompson       if (use_3d_slices) {
449*4b3e95d5SJeremy L Thompson         // Accumulator for gradient slices
450*4b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n";
451*4b3e95d5SJeremy L Thompson         code << "    for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n";
452*4b3e95d5SJeremy L Thompson         code << "      r_q" << var_suffix << "[i] = 0.0;\n";
453*4b3e95d5SJeremy L Thompson         code << "    }\n";
454*4b3e95d5SJeremy L Thompson       } else {
455*4b3e95d5SJeremy L Thompson         code << "    CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n";
456*4b3e95d5SJeremy L Thompson       }
457*4b3e95d5SJeremy L Thompson     }
458*4b3e95d5SJeremy L Thompson   }
459*4b3e95d5SJeremy L Thompson 
460*4b3e95d5SJeremy L Thompson   // We treat quadrature points per slice in 3d to save registers
461*4b3e95d5SJeremy L Thompson   if (use_3d_slices) {
462*4b3e95d5SJeremy L Thompson     code << "\n    // Note: Using planes of 3D elements\n";
463*4b3e95d5SJeremy L Thompson     code << "#pragma unroll\n";
464*4b3e95d5SJeremy L Thompson     code << "    for (CeedInt q = 0; q < " << Q_name << "; q++) {\n";
465*4b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
466*4b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
467*4b3e95d5SJeremy L Thompson       std::string var_suffix = "_in_" + std::to_string(i);
468*4b3e95d5SJeremy L Thompson 
469*4b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
470*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
471*4b3e95d5SJeremy L Thompson       // Basis action
472*4b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
473*4b3e95d5SJeremy L Thompson       switch (eval_mode) {
474*4b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
475*4b3e95d5SJeremy L Thompson           bool is_strided;
476*4b3e95d5SJeremy L Thompson 
477*4b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
478*4b3e95d5SJeremy L Thompson 
479*4b3e95d5SJeremy L Thompson           CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr));
480*4b3e95d5SJeremy L Thompson           CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided));
481*4b3e95d5SJeremy L Thompson           if (is_strided) {
482*4b3e95d5SJeremy L Thompson             bool    has_backend_strides;
483*4b3e95d5SJeremy L Thompson             CeedInt num_elem, elem_size;
484*4b3e95d5SJeremy L Thompson 
485*4b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size));
486*4b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides));
487*4b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem));
488*4b3e95d5SJeremy L Thompson             CeedInt strides[3] = {1, elem_size * num_elem, elem_size};
489*4b3e95d5SJeremy L Thompson 
490*4b3e95d5SJeremy L Thompson             if (!has_backend_strides) {
491*4b3e95d5SJeremy L Thompson               CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides));
492*4b3e95d5SJeremy L Thompson             }
493*4b3e95d5SJeremy L Thompson             code << "      // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n";
494*4b3e95d5SJeremy L Thompson             code << "      readSliceQuadsStrided3d<num_comp" << var_suffix << ", " << Q_name << "," << strides[0] << "," << strides[1] << ","
495*4b3e95d5SJeremy L Thompson                  << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n";
496*4b3e95d5SJeremy L Thompson           } else {
497*4b3e95d5SJeremy L Thompson             CeedSize                 l_size = 0;
498*4b3e95d5SJeremy L Thompson             CeedInt                  comp_stride;
499*4b3e95d5SJeremy L Thompson             CeedElemRestriction_Hip *rstr_data;
500*4b3e95d5SJeremy L Thompson 
501*4b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size));
502*4b3e95d5SJeremy L Thompson             code << "      const CeedInt l_size" << var_suffix << " = " << l_size << ";\n";
503*4b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride));
504*4b3e95d5SJeremy L Thompson             code << "      // CompStride: " << comp_stride << "\n";
505*4b3e95d5SJeremy L Thompson             CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data));
506*4b3e95d5SJeremy L Thompson             data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets;
507*4b3e95d5SJeremy L Thompson             code << "      readSliceQuadsOffset3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix
508*4b3e95d5SJeremy L Thompson                  << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n";
509*4b3e95d5SJeremy L Thompson           }
510*4b3e95d5SJeremy L Thompson           break;
511*4b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
512*4b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
513*4b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n";
514*4b3e95d5SJeremy L Thompson           code << "        r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n";
515*4b3e95d5SJeremy L Thompson           code << "      }\n";
516*4b3e95d5SJeremy L Thompson           break;
517*4b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
518*4b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
519*4b3e95d5SJeremy L Thompson           code << "      gradCollo3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix << ", r_s"
520*4b3e95d5SJeremy L Thompson                << var_suffix << ");\n";
521*4b3e95d5SJeremy L Thompson           break;
522*4b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
523*4b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[1];\n";
524*4b3e95d5SJeremy L Thompson           code << "      r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n";
525*4b3e95d5SJeremy L Thompson           break;  // No action
526*4b3e95d5SJeremy L Thompson                   // LCOV_EXCL_START
527*4b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
528*4b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
529*4b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
530*4b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
531*4b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
532*4b3e95d5SJeremy L Thompson       }
533*4b3e95d5SJeremy L Thompson     }
534*4b3e95d5SJeremy L Thompson     code << "\n      // -- Output fields\n";
535*4b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
536*4b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
537*4b3e95d5SJeremy L Thompson 
538*4b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
539*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
540*4b3e95d5SJeremy L Thompson       // Basis action
541*4b3e95d5SJeremy L Thompson       switch (eval_mode) {
542*4b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
543*4b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
544*4b3e95d5SJeremy L Thompson           break;  // No action
545*4b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
546*4b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n";
547*4b3e95d5SJeremy L Thompson           break;
548*4b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
549*4b3e95d5SJeremy L Thompson           code << "      CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n";
550*4b3e95d5SJeremy L Thompson           break;
551*4b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
552*4b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
553*4b3e95d5SJeremy L Thompson           break;  // Should not occur
554*4b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
555*4b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
556*4b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
557*4b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
558*4b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
559*4b3e95d5SJeremy L Thompson       }
560*4b3e95d5SJeremy L Thompson     }
561*4b3e95d5SJeremy L Thompson   } else {
562*4b3e95d5SJeremy L Thompson     code << "\n    // Note: Using full elements\n";
563*4b3e95d5SJeremy L Thompson     code << "    {\n";
564*4b3e95d5SJeremy L Thompson     code << "      // -- Input fields\n";
565*4b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_input_fields; i++) {
566*4b3e95d5SJeremy L Thompson       code << "      // ---- Input field " << i << "\n";
567*4b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n";
568*4b3e95d5SJeremy L Thompson     }
569*4b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
570*4b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
571*4b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
572*4b3e95d5SJeremy L Thompson       code << "      CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n";
573*4b3e95d5SJeremy L Thompson     }
574*4b3e95d5SJeremy L Thompson   }
575*4b3e95d5SJeremy L Thompson 
576*4b3e95d5SJeremy L Thompson   // Input and output buffers
577*4b3e95d5SJeremy L Thompson   code << "\n      // -- QFunction inputs and outputs\n";
578*4b3e95d5SJeremy L Thompson   code << "      // ---- Inputs\n";
579*4b3e95d5SJeremy L Thompson   code << "      CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
580*4b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
581*4b3e95d5SJeremy L Thompson     code << "      // ------ Input field " << i << "\n";
582*4b3e95d5SJeremy L Thompson     code << "      inputs[" << i << "] = r_s_in_" << i << ";\n";
583*4b3e95d5SJeremy L Thompson   }
584*4b3e95d5SJeremy L Thompson   code << "      // ---- Outputs\n";
585*4b3e95d5SJeremy L Thompson   code << "      CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
586*4b3e95d5SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
587*4b3e95d5SJeremy L Thompson     code << "      // ------ Output field " << i << "\n";
588*4b3e95d5SJeremy L Thompson     code << "      outputs[" << i << "] = r_s_out_" << i << ";\n";
589*4b3e95d5SJeremy L Thompson   }
590*4b3e95d5SJeremy L Thompson 
591*4b3e95d5SJeremy L Thompson   // Apply QFunction
592*4b3e95d5SJeremy L Thompson   code << "\n      // -- Apply QFunction\n";
593*4b3e95d5SJeremy L Thompson   code << "      " << qfunction_name << "(ctx, ";
594*4b3e95d5SJeremy L Thompson   if (dim != 3 || use_3d_slices) {
595*4b3e95d5SJeremy L Thompson     code << "1";
596*4b3e95d5SJeremy L Thompson   } else {
597*4b3e95d5SJeremy L Thompson     code << "Q_1d";
598*4b3e95d5SJeremy L Thompson   }
599*4b3e95d5SJeremy L Thompson   code << ", inputs, outputs);\n";
600*4b3e95d5SJeremy L Thompson 
601*4b3e95d5SJeremy L Thompson   // Copy or apply transpose grad, if needed
602*4b3e95d5SJeremy L Thompson   if (use_3d_slices) {
603*4b3e95d5SJeremy L Thompson     code << "      // -- Output fields\n";
604*4b3e95d5SJeremy L Thompson     for (CeedInt i = 0; i < num_output_fields; i++) {
605*4b3e95d5SJeremy L Thompson       std::string var_suffix = "_out_" + std::to_string(i);
606*4b3e95d5SJeremy L Thompson       std::string P_name     = "P_1d" + var_suffix;
607*4b3e95d5SJeremy L Thompson 
608*4b3e95d5SJeremy L Thompson       code << "      // ---- Output field " << i << "\n";
609*4b3e95d5SJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
610*4b3e95d5SJeremy L Thompson       // Basis action
611*4b3e95d5SJeremy L Thompson       code << "      // EvalMode: " << CeedEvalModes[eval_mode] << "\n";
612*4b3e95d5SJeremy L Thompson       switch (eval_mode) {
613*4b3e95d5SJeremy L Thompson         case CEED_EVAL_NONE:
614*4b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
615*4b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
616*4b3e95d5SJeremy L Thompson           code << "      }\n";
617*4b3e95d5SJeremy L Thompson           break;  // No action
618*4b3e95d5SJeremy L Thompson         case CEED_EVAL_INTERP:
619*4b3e95d5SJeremy L Thompson           code << "      for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n";
620*4b3e95d5SJeremy L Thompson           code << "        r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n";
621*4b3e95d5SJeremy L Thompson           code << "      }\n";
622*4b3e95d5SJeremy L Thompson           break;
623*4b3e95d5SJeremy L Thompson         case CEED_EVAL_GRAD:
624*4b3e95d5SJeremy L Thompson           code << "      gradColloTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G" << var_suffix
625*4b3e95d5SJeremy L Thompson                << ", r_q" << var_suffix << ");\n";
626*4b3e95d5SJeremy L Thompson           break;
627*4b3e95d5SJeremy L Thompson           // LCOV_EXCL_START
628*4b3e95d5SJeremy L Thompson         case CEED_EVAL_WEIGHT:
629*4b3e95d5SJeremy L Thompson           break;  // Should not occur
630*4b3e95d5SJeremy L Thompson         case CEED_EVAL_DIV:
631*4b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
632*4b3e95d5SJeremy L Thompson         case CEED_EVAL_CURL:
633*4b3e95d5SJeremy L Thompson           break;  // TODO: Not implemented
634*4b3e95d5SJeremy L Thompson                   // LCOV_EXCL_STOP
635*4b3e95d5SJeremy L Thompson       }
636*4b3e95d5SJeremy L Thompson     }
637*4b3e95d5SJeremy L Thompson   }
638*4b3e95d5SJeremy L Thompson   code << "    }\n";
639*4b3e95d5SJeremy L Thompson   return CEED_ERROR_SUCCESS;
640*4b3e95d5SJeremy L Thompson }
641*4b3e95d5SJeremy L Thompson 
642*4b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------
6439e201c85SYohann // Build single operator kernel
6447d8d0e25Snbeams //------------------------------------------------------------------------------
645eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op) {
646*4b3e95d5SJeremy L Thompson   bool                   is_tensor = true, use_3d_slices = false;
6477d8d0e25Snbeams   Ceed                   ceed;
648*4b3e95d5SJeremy L Thompson   CeedInt                Q_1d, num_input_fields, num_output_fields, dim = 1;
649b7453713SJeremy L Thompson   CeedQFunctionField    *qf_input_fields, *qf_output_fields;
650b7453713SJeremy L Thompson   CeedQFunction_Hip_gen *qf_data;
651b7453713SJeremy L Thompson   CeedQFunction          qf;
652b7453713SJeremy L Thompson   CeedOperatorField     *op_input_fields, *op_output_fields;
653b7453713SJeremy L Thompson   CeedOperator_Hip_gen  *data;
654*4b3e95d5SJeremy L Thompson   std::ostringstream     code;
655*4b3e95d5SJeremy L Thompson 
656*4b3e95d5SJeremy L Thompson   {
657*4b3e95d5SJeremy L Thompson     bool is_setup_done;
658b7453713SJeremy L Thompson 
659b7453713SJeremy L Thompson     CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done));
660b7453713SJeremy L Thompson     if (is_setup_done) return CEED_ERROR_SUCCESS;
661*4b3e95d5SJeremy L Thompson   }
662b7453713SJeremy L Thompson 
663b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
664b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
665b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
666b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
667b7453713SJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
668b7453713SJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
6697d8d0e25Snbeams 
670*4b3e95d5SJeremy L Thompson   // Get operator data
671*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelData_Hip_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields,
672*4b3e95d5SJeremy L Thompson                                                       qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices));
673*4b3e95d5SJeremy L Thompson   if (dim == 0) dim = 1;
674*4b3e95d5SJeremy L Thompson   data->dim = dim;
675*4b3e95d5SJeremy L Thompson   if (Q_1d == 0) {
676*4b3e95d5SJeremy L Thompson     CeedInt Q;
677*4b3e95d5SJeremy L Thompson 
678*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q));
679*4b3e95d5SJeremy L Thompson     Q_1d = Q;
680*4b3e95d5SJeremy L Thompson   }
681*4b3e95d5SJeremy L Thompson   data->Q_1d = Q_1d;
682*4b3e95d5SJeremy L Thompson 
6830b454692Sjeremylt   // Check for restriction only identity operator
684*4b3e95d5SJeremy L Thompson   {
685*4b3e95d5SJeremy L Thompson     bool is_identity_qf;
686*4b3e95d5SJeremy L Thompson 
6872b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf));
6880b454692Sjeremylt     if (is_identity_qf) {
6899e201c85SYohann       CeedEvalMode eval_mode_in, eval_mode_out;
690b7453713SJeremy L Thompson 
6912b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in));
6922b730f8bSJeremy L Thompson       CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out));
6936574a04fSJeremy L Thompson       CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND,
6946574a04fSJeremy L Thompson                 "Backend does not implement restriction only identity operators");
6950b454692Sjeremylt     }
696*4b3e95d5SJeremy L Thompson   }
697b2165e7aSSebastian Grimberg 
698b2165e7aSSebastian Grimberg   // Load basis source files
699*4b3e95d5SJeremy L Thompson   // TODO: Add non-tensor, AtPoints
7009e201c85SYohann   {
70122070f95SJeremy L Thompson     char       *tensor_basis_kernel_source;
70222070f95SJeremy L Thompson     const char *tensor_basis_kernel_path;
703b7453713SJeremy L Thompson 
7042b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-shared-basis-tensor-templates.h", &tensor_basis_kernel_path));
70523d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Tensor Basis Kernel Source -----\n");
7062b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source));
7079e201c85SYohann     code << tensor_basis_kernel_source;
7082b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_path));
7092b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&tensor_basis_kernel_source));
7109e201c85SYohann   }
7119e201c85SYohann   {
71222070f95SJeremy L Thompson     char       *hip_gen_template_source;
71322070f95SJeremy L Thompson     const char *hip_gen_template_path;
714b7453713SJeremy L Thompson 
7152b730f8bSJeremy L Thompson     CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-gen-templates.h", &hip_gen_template_path));
71623d4529eSJeremy L Thompson     CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Hip-Gen Template Source -----\n");
7172b730f8bSJeremy L Thompson     CeedCallBackend(CeedLoadSourceToBuffer(ceed, hip_gen_template_path, &hip_gen_template_source));
7189e201c85SYohann     code << hip_gen_template_source;
7192b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&hip_gen_template_path));
7202b730f8bSJeremy L Thompson     CeedCallBackend(CeedFree(&hip_gen_template_source));
7219e201c85SYohann   }
7227d8d0e25Snbeams 
723*4b3e95d5SJeremy L Thompson   // Get QFunction name
724*4b3e95d5SJeremy L Thompson   std::string qfunction_name(qf_data->qfunction_name);
725*4b3e95d5SJeremy L Thompson   std::string operator_name;
726*4b3e95d5SJeremy L Thompson 
72709095acaSJeremy L Thompson   operator_name = "CeedKernelHipGenOperator_" + qfunction_name;
7287d8d0e25Snbeams 
7299e201c85SYohann   // Define CEED_Q_VLA
7309e201c85SYohann   code << "\n#undef CEED_Q_VLA\n";
731*4b3e95d5SJeremy L Thompson   if (dim != 3 || use_3d_slices) {
7329e201c85SYohann     code << "#define CEED_Q_VLA 1\n\n";
7339e201c85SYohann   } else {
7349e201c85SYohann     code << "#define CEED_Q_VLA " << Q_1d << "\n\n";
7359e201c85SYohann   }
7369e201c85SYohann 
737*4b3e95d5SJeremy L Thompson   // Add user QFunction source
738*4b3e95d5SJeremy L Thompson   {
739*4b3e95d5SJeremy L Thompson     std::string qfunction_source(qf_data->qfunction_source);
740*4b3e95d5SJeremy L Thompson 
74109095acaSJeremy L Thompson     code << qfunction_source;
742*4b3e95d5SJeremy L Thompson   }
7437d8d0e25Snbeams 
7447d8d0e25Snbeams   // Setup
7457d8d0e25Snbeams   code << "\n// -----------------------------------------------------------------------------\n";
746*4b3e95d5SJeremy L Thompson   code << "// Operator Kernel\n";
747*4b3e95d5SJeremy L Thompson   code << "// \n";
748*4b3e95d5SJeremy L Thompson   code << "// d_[in,out]_i:   CeedVector device array\n";
749*4b3e95d5SJeremy L Thompson   code << "// r_[in,out]_e_i: Element vector register\n";
750*4b3e95d5SJeremy L Thompson   code << "// r_[in,out]_q_i: Quadrature space vector register\n";
751*4b3e95d5SJeremy L Thompson   code << "// r_[in,out]_s_i: Quadrature space slice  vector register\n";
752*4b3e95d5SJeremy L Thompson   code << "// \n";
753*4b3e95d5SJeremy L Thompson   code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n";
754*4b3e95d5SJeremy L Thompson   code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n";
755*4b3e95d5SJeremy L Thompson   code << "// -----------------------------------------------------------------------------\n";
756b3e1519bSnbeams   code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
7572b730f8bSJeremy L Thompson   code << "__global__ void " << operator_name
7582b730f8bSJeremy L Thompson        << "(CeedInt num_elem, void* ctx, FieldsInt_Hip indices, Fields_Hip fields, Fields_Hip B, Fields_Hip G, CeedScalar* W) {\n";
759*4b3e95d5SJeremy L Thompson 
760*4b3e95d5SJeremy L Thompson   // Scratch buffers
7619e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
762*4b3e95d5SJeremy L Thompson     CeedEvalMode eval_mode;
763*4b3e95d5SJeremy L Thompson 
7642b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
7659e201c85SYohann     if (eval_mode != CEED_EVAL_WEIGHT) {  // Skip CEED_EVAL_WEIGHT
766*4b3e95d5SJeremy L Thompson       code << "  const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n";
7677d8d0e25Snbeams     }
7687d8d0e25Snbeams   }
7699e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
770*4b3e95d5SJeremy L Thompson     code << "  CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n";
7717d8d0e25Snbeams   }
7727d8d0e25Snbeams 
7739e201c85SYohann   code << "  const CeedInt dim = " << dim << ";\n";
7749e201c85SYohann   code << "  const CeedInt Q_1d = " << Q_1d << ";\n";
7757d8d0e25Snbeams 
776*4b3e95d5SJeremy L Thompson   // Shared data
777*4b3e95d5SJeremy L Thompson   code << "  extern __shared__ CeedScalar slice[];\n";
7789e201c85SYohann   code << "  SharedData_Hip data;\n";
7799e201c85SYohann   code << "  data.t_id_x = threadIdx.x;\n";
7809e201c85SYohann   code << "  data.t_id_y = threadIdx.y;\n";
7819e201c85SYohann   code << "  data.t_id_z = threadIdx.z;\n";
7829e201c85SYohann   code << "  data.t_id  = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n";
7839e201c85SYohann   code << "  data.slice = slice + data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n";
7847d8d0e25Snbeams 
7857d8d0e25Snbeams   // Initialize constants, and matrices B and G
786*4b3e95d5SJeremy L Thompson   code << "\n  // Input field constants and basis data\n";
7879e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
788*4b3e95d5SJeremy L Thompson     CeedCall(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_input_fields[i], qf_input_fields[i], Q_1d, true, use_3d_slices));
7897d8d0e25Snbeams   }
790*4b3e95d5SJeremy L Thompson   code << "\n  // Output field constants and basis data\n";
7919e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
792*4b3e95d5SJeremy L Thompson     CeedCall(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices));
793*4b3e95d5SJeremy L Thompson   }
7947d8d0e25Snbeams 
795*4b3e95d5SJeremy L Thompson   // Loop over all elements
796*4b3e95d5SJeremy L Thompson   code << "\n  // Element loop\n";
7977d8d0e25Snbeams   code << "  __syncthreads();\n";
7989e201c85SYohann   code << "  for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n";
799*4b3e95d5SJeremy L Thompson 
800*4b3e95d5SJeremy L Thompson   // -- Input restriction and basis
801*4b3e95d5SJeremy L Thompson   code << "    // -- Input field restrictions and basis actions\n";
8029e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
803*4b3e95d5SJeremy L Thompson     code << "    // ---- Input field " << i << "\n";
8047d8d0e25Snbeams 
805*4b3e95d5SJeremy L Thompson     // ---- Restriction
806*4b3e95d5SJeremy L Thompson     CeedCallBackend(
807*4b3e95d5SJeremy L Thompson         CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, dim, op_input_fields[i], qf_input_fields[i], Q_1d, true, use_3d_slices));
808b7453713SJeremy L Thompson 
809*4b3e95d5SJeremy L Thompson     // ---- Basis action
810*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, dim, op_input_fields[i], qf_input_fields[i], Q_1d, true, use_3d_slices));
8117d8d0e25Snbeams   }
8127d8d0e25Snbeams 
813*4b3e95d5SJeremy L Thompson   // -- Q function
814*4b3e95d5SJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernelQFunction_Hip_gen(code, data, dim, num_input_fields, op_input_fields, qf_input_fields, num_output_fields,
815*4b3e95d5SJeremy L Thompson                                                            op_output_fields, qf_output_fields, qfunction_name, Q_1d, use_3d_slices));
8167d8d0e25Snbeams 
817*4b3e95d5SJeremy L Thompson   // -- Output basis and restriction
818*4b3e95d5SJeremy L Thompson   code << "\n    // -- Output field basis action and restrictions\n";
8199e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
820*4b3e95d5SJeremy L Thompson     code << "    // ---- Output field " << i << "\n";
821b7453713SJeremy L Thompson 
822*4b3e95d5SJeremy L Thompson     // ---- Basis action
823*4b3e95d5SJeremy L Thompson     CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices));
8247d8d0e25Snbeams 
825*4b3e95d5SJeremy L Thompson     // ---- Restriction
826*4b3e95d5SJeremy L Thompson     CeedCallBackend(
827*4b3e95d5SJeremy L Thompson         CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices));
8287d8d0e25Snbeams   }
8297d8d0e25Snbeams 
830*4b3e95d5SJeremy L Thompson   // Close loop and function
8317d8d0e25Snbeams   code << "  }\n";
8327d8d0e25Snbeams   code << "}\n";
8337d8d0e25Snbeams   code << "// -----------------------------------------------------------------------------\n\n";
8347d8d0e25Snbeams 
8357d8d0e25Snbeams   // View kernel for debugging
83623d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated Operator Kernels:\n");
8373f21f6b1SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
8387d8d0e25Snbeams 
839539ec17dSJeremy L Thompson   CeedInt block_sizes[3] = {0, 0, 0};
8409e201c85SYohann   CeedInt num_elem;
841b7453713SJeremy L Thompson 
8422b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
8432b730f8bSJeremy L Thompson   CeedCallBackend(BlockGridCalculate_Hip_gen(dim, num_elem, data->max_P_1d, Q_1d, block_sizes));
844eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE",
8452b730f8bSJeremy L Thompson                                   block_sizes[0] * block_sizes[1] * block_sizes[2]));
846eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op));
8477d8d0e25Snbeams 
8482b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetSetupDone(op));
849e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
8507d8d0e25Snbeams }
8512a86cc9dSSebastian Grimberg 
8527d8d0e25Snbeams //------------------------------------------------------------------------------
853