1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 #define CEED_DEBUG_COLOR 12 9 10 #include <ceed.h> 11 #include <ceed/backend.h> 12 #include <ceed/jit-tools.h> 13 14 #include <iostream> 15 #include <sstream> 16 #include <string> 17 18 #include "../hip-ref/ceed-hip-ref.h" 19 #include "../hip-shared/ceed-hip-shared.h" 20 #include "../hip/ceed-hip-common.h" 21 #include "../hip/ceed-hip-compile.h" 22 #include "ceed-hip-gen.h" 23 24 //------------------------------------------------------------------------------ 25 // Calculate the block size used for launching the operator kernel 26 //------------------------------------------------------------------------------ 27 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 const CeedInt thread1d = CeedIntMax(Q_1d, P_1d); 29 if (dim == 1) { 30 CeedInt elems_per_block = 64 * thread1d > 256 ? 256 / thread1d : 64; 31 32 elems_per_block = elems_per_block > 0 ? elems_per_block : 1; 33 block_sizes[0] = thread1d; 34 block_sizes[1] = 1; 35 block_sizes[2] = elems_per_block; 36 } else if (dim == 2) { 37 const CeedInt elems_per_block = thread1d < 4 ? 16 : 2; 38 39 block_sizes[0] = thread1d; 40 block_sizes[1] = thread1d; 41 block_sizes[2] = elems_per_block; 42 } else if (dim == 3) { 43 const CeedInt elems_per_block = thread1d < 6 ? 4 : (thread1d < 8 ? 2 : 1); 44 45 block_sizes[0] = thread1d; 46 block_sizes[1] = thread1d; 47 block_sizes[2] = elems_per_block; 48 } 49 return CEED_ERROR_SUCCESS; 50 } 51 52 //------------------------------------------------------------------------------ 53 // Determine type of operator 54 //------------------------------------------------------------------------------ 55 static int CeedOperatorBuildKernelData_Hip_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields, 56 CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields, 57 CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor, 58 bool *use_3d_slices) { 59 // Find dim, P_1d, Q_1d 60 *max_P_1d = 0; 61 *Q_1d = 0; 62 *dim = 0; 63 *is_tensor = true; 64 for (CeedInt i = 0; i < num_input_fields; i++) { 65 CeedBasis basis; 66 67 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 68 if (basis != CEED_BASIS_NONE) { 69 bool is_field_tensor; 70 CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0; 71 72 // Collect dim, P_1d, and Q_1d 73 CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 74 CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 75 *is_tensor = *is_tensor && is_field_tensor; 76 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 77 *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 78 CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 79 CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 80 *dim = field_dim; 81 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 82 CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 83 *Q_1d = field_Q_1d; 84 } 85 } 86 for (CeedInt i = 0; i < num_output_fields; i++) { 87 CeedBasis basis; 88 89 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 90 if (basis != CEED_BASIS_NONE) { 91 bool is_field_tensor; 92 CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0; 93 94 // Collect dim, P_1d, and Q_1d 95 CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 96 CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 97 *is_tensor = *is_tensor && is_field_tensor; 98 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 99 *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 100 CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 101 CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 102 *dim = field_dim; 103 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 104 CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 105 *Q_1d = field_Q_1d; 106 } 107 } 108 109 // Only use 3D collocated gradient parallelization strategy when gradient is computed 110 *use_3d_slices = false; 111 if (*dim == 3) { 112 bool was_grad_found = false; 113 114 for (CeedInt i = 0; i < num_input_fields; i++) { 115 CeedEvalMode eval_mode; 116 117 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 118 if (eval_mode == CEED_EVAL_GRAD) { 119 CeedBasis_Hip_shared *basis_data; 120 CeedBasis basis; 121 122 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 123 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 124 *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 125 was_grad_found = true; 126 } 127 } 128 for (CeedInt i = 0; i < num_output_fields; i++) { 129 CeedEvalMode eval_mode; 130 131 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 132 if (eval_mode == CEED_EVAL_GRAD) { 133 CeedBasis_Hip_shared *basis_data; 134 CeedBasis basis; 135 136 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 137 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 138 *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 139 was_grad_found = true; 140 } 141 } 142 } 143 return CEED_ERROR_SUCCESS; 144 } 145 146 //------------------------------------------------------------------------------ 147 // Setup fields 148 //------------------------------------------------------------------------------ 149 static int CeedOperatorBuildKernelFieldData_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedOperatorField op_field, 150 CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool use_3d_slices) { 151 std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 152 std::string P_name = "P_1d" + var_suffix, Q_name = "Q_1d"; 153 std::string option_name = (is_input ? "inputs" : "outputs"); 154 CeedEvalMode eval_mode = CEED_EVAL_NONE; 155 CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 156 CeedElemRestriction elem_rstr; 157 CeedBasis_Hip_shared *basis_data; 158 CeedBasis basis; 159 160 code << " // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n"; 161 162 // Get field data 163 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 164 if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 165 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 166 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 167 } 168 CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 169 if (basis != CEED_BASIS_NONE) { 170 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 171 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 172 } 173 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 174 175 // Set field constants 176 if (eval_mode != CEED_EVAL_WEIGHT) { 177 code << " const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n"; 178 code << " const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n"; 179 } 180 181 // Load basis data 182 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 183 switch (eval_mode) { 184 case CEED_EVAL_NONE: 185 break; 186 case CEED_EVAL_INTERP: 187 if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 188 else data->B.outputs[i] = basis_data->d_interp_1d; 189 code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n"; 190 code << " loadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 191 break; 192 case CEED_EVAL_GRAD: 193 if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 194 else data->B.outputs[i] = basis_data->d_interp_1d; 195 code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n"; 196 code << " loadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 197 if (use_3d_slices) { 198 if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d; 199 else data->G.outputs[i] = basis_data->d_collo_grad_1d; 200 code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n"; 201 code << " loadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 202 } else { 203 bool has_collo_grad = basis_data->d_collo_grad_1d; 204 205 if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 206 else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 207 if (has_collo_grad) { 208 code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n"; 209 code << " loadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 210 } else { 211 code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * P_1d << "];\n"; 212 code << " loadMatrix<" << P_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 213 } 214 } 215 break; 216 case CEED_EVAL_WEIGHT: 217 break; // No action 218 // LCOV_EXCL_START 219 case CEED_EVAL_DIV: 220 break; // TODO: Not implemented 221 case CEED_EVAL_CURL: 222 break; // TODO: Not implemented 223 // LCOV_EXCL_STOP 224 } 225 return CEED_ERROR_SUCCESS; 226 } 227 228 //------------------------------------------------------------------------------ 229 // Restriction 230 //------------------------------------------------------------------------------ 231 static int CeedOperatorBuildKernelRestriction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim, 232 CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field, 233 CeedInt Q_1d, bool is_input, bool use_3d_slices) { 234 std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 235 std::string P_name = "P_1d" + var_suffix; 236 CeedEvalMode eval_mode = CEED_EVAL_NONE; 237 CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 238 CeedSize l_size; 239 CeedElemRestriction_Hip *rstr_data; 240 CeedElemRestriction elem_rstr; 241 CeedBasis basis; 242 243 // Get field data 244 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 245 if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 246 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 247 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 248 CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 249 } 250 CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 251 if (basis != CEED_BASIS_NONE) { 252 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 253 } 254 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 255 256 // Restriction 257 if (is_input) { 258 // Input 259 // Input 260 if (field_input_buffer[i] != i) { 261 std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]); 262 263 // Restriction was already done for previous input 264 code << " CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n"; 265 } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices)) { 266 bool is_strided; 267 268 if (eval_mode == CEED_EVAL_NONE) { 269 // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated 270 code << " CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n"; 271 } else { 272 // Otherwise we're using the scratch space 273 code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 274 } 275 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 276 if (!is_strided) { 277 CeedInt comp_stride; 278 279 CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 280 code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 281 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 282 code << " // CompStride: " << comp_stride << "\n"; 283 data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 284 code << " readDofsOffset" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" << var_suffix 285 << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n"; 286 } else { 287 bool has_backend_strides; 288 CeedInt num_elem; 289 290 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 291 CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 292 CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 293 294 if (!has_backend_strides) { 295 CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 296 } 297 code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 298 code << " readDofsStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << "," 299 << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n"; 300 } 301 } 302 } else { 303 // Output 304 bool is_strided; 305 306 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 307 if (!is_strided) { 308 CeedInt comp_stride; 309 310 CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 311 code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 312 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 313 code << " // CompStride: " << comp_stride << "\n"; 314 data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets; 315 code << " writeDofsOffset" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" << var_suffix 316 << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n"; 317 } else { 318 bool has_backend_strides; 319 CeedInt num_elem; 320 321 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 322 CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 323 CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 324 325 if (!has_backend_strides) { 326 CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 327 } 328 code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 329 code << " writeDofsStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << "," 330 << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n"; 331 } 332 } 333 return CEED_ERROR_SUCCESS; 334 } 335 336 //------------------------------------------------------------------------------ 337 // Basis 338 //------------------------------------------------------------------------------ 339 static int CeedOperatorBuildKernelBasis_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim, 340 CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, 341 bool use_3d_slices) { 342 std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 343 std::string P_name = "P_1d" + var_suffix, Q_name = "Q_1d"; 344 CeedEvalMode eval_mode = CEED_EVAL_NONE; 345 CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 346 CeedElemRestriction elem_rstr; 347 CeedBasis basis; 348 349 // Get field data 350 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 351 if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 352 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 353 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 354 } 355 CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 356 if (basis != CEED_BASIS_NONE) { 357 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 358 } 359 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 360 361 // Basis 362 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 363 if (is_input) { 364 switch (eval_mode) { 365 case CEED_EVAL_NONE: 366 if (!use_3d_slices) { 367 code << " CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n"; 368 } 369 break; 370 case CEED_EVAL_INTERP: 371 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 372 code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name 373 << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n"; 374 break; 375 case CEED_EVAL_GRAD: 376 if (use_3d_slices) { 377 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 378 code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name 379 << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n"; 380 } else { 381 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n"; 382 code << " Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" << var_suffix 383 << ", P_1d" << var_suffix << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q" 384 << var_suffix << ");\n"; 385 } 386 break; 387 case CEED_EVAL_WEIGHT: { 388 CeedBasis_Hip_shared *basis_data; 389 390 code << " CeedScalar r_q" << var_suffix << "[" << Q_name << "];\n"; 391 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 392 data->W = basis_data->d_q_weight_1d; 393 code << " Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n"; 394 break; 395 } 396 // LCOV_EXCL_START 397 case CEED_EVAL_DIV: 398 case CEED_EVAL_CURL: 399 break; // TODO: Not implemented 400 // LCOV_EXCL_STOP 401 } 402 } else { 403 switch (eval_mode) { 404 case CEED_EVAL_NONE: 405 code << " CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n"; 406 break; // No action 407 case CEED_EVAL_INTERP: 408 code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 409 code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name 410 << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 411 break; 412 case CEED_EVAL_GRAD: 413 code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 414 if (use_3d_slices) { 415 code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name 416 << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 417 } else { 418 code << " GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" 419 << var_suffix << ", " << P_name << "," << Q_name << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix 420 << ", r_e" << var_suffix << ");\n"; 421 } 422 break; 423 // LCOV_EXCL_START 424 case CEED_EVAL_WEIGHT: 425 break; // Should not occur 426 case CEED_EVAL_DIV: 427 case CEED_EVAL_CURL: 428 break; // TODO: Not implemented 429 // LCOV_EXCL_STOP 430 } 431 } 432 return CEED_ERROR_SUCCESS; 433 } 434 435 //------------------------------------------------------------------------------ 436 // QFunction 437 //------------------------------------------------------------------------------ 438 static int CeedOperatorBuildKernelQFunction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt dim, CeedInt num_input_fields, 439 CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields, 440 CeedInt num_output_fields, CeedOperatorField *op_output_fields, 441 CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d, 442 bool use_3d_slices) { 443 std::string Q_name = "Q_1d"; 444 CeedEvalMode eval_mode = CEED_EVAL_NONE; 445 CeedElemRestriction elem_rstr; 446 447 // Setup output arays 448 code << "\n // -- Output field setup\n"; 449 for (CeedInt i = 0; i < num_output_fields; i++) { 450 std::string var_suffix = "_out_" + std::to_string(i); 451 452 code << " // ---- Output field " << i << "\n"; 453 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 454 if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) { 455 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 456 } 457 if (eval_mode == CEED_EVAL_GRAD) { 458 if (use_3d_slices) { 459 // Accumulator for gradient slices 460 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 461 code << " for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n"; 462 code << " r_q" << var_suffix << "[i] = 0.0;\n"; 463 code << " }\n"; 464 } else { 465 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n"; 466 } 467 } 468 } 469 470 // We treat quadrature points per slice in 3d to save registers 471 if (use_3d_slices) { 472 code << "\n // Note: Using planes of 3D elements\n"; 473 code << " #pragma unroll\n"; 474 code << " for (CeedInt q = 0; q < " << Q_name << "; q++) {\n"; 475 code << " // -- Input fields\n"; 476 for (CeedInt i = 0; i < num_input_fields; i++) { 477 std::string var_suffix = "_in_" + std::to_string(i); 478 479 code << " // ---- Input field " << i << "\n"; 480 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 481 // Basis action 482 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 483 switch (eval_mode) { 484 case CEED_EVAL_NONE: 485 bool is_strided; 486 487 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 488 489 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 490 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 491 if (is_strided) { 492 bool has_backend_strides; 493 CeedInt num_elem, elem_size; 494 495 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 496 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 497 CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 498 CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 499 500 if (!has_backend_strides) { 501 CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 502 } 503 code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 504 code << " readSliceQuadsStrided3d<num_comp" << var_suffix << ", " << Q_name << "," << strides[0] << "," << strides[1] << "," 505 << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n"; 506 } else { 507 CeedSize l_size = 0; 508 CeedInt comp_stride; 509 CeedElemRestriction_Hip *rstr_data; 510 511 CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 512 code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 513 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 514 code << " // CompStride: " << comp_stride << "\n"; 515 CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 516 data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 517 code << " readSliceQuadsOffset3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix 518 << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n"; 519 } 520 break; 521 case CEED_EVAL_INTERP: 522 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 523 code << " for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n"; 524 code << " r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n"; 525 code << " }\n"; 526 break; 527 case CEED_EVAL_GRAD: 528 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 529 code << " gradCollo3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix << ", r_s" 530 << var_suffix << ");\n"; 531 break; 532 case CEED_EVAL_WEIGHT: 533 code << " CeedScalar r_s" << var_suffix << "[1];\n"; 534 code << " r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n"; 535 break; // No action 536 // LCOV_EXCL_START 537 case CEED_EVAL_DIV: 538 break; // TODO: Not implemented 539 case CEED_EVAL_CURL: 540 break; // TODO: Not implemented 541 // LCOV_EXCL_STOP 542 } 543 } 544 code << "\n // -- Output fields\n"; 545 for (CeedInt i = 0; i < num_output_fields; i++) { 546 std::string var_suffix = "_out_" + std::to_string(i); 547 548 code << " // ---- Output field " << i << "\n"; 549 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 550 // Basis action 551 switch (eval_mode) { 552 case CEED_EVAL_NONE: 553 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 554 break; // No action 555 case CEED_EVAL_INTERP: 556 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 557 break; 558 case CEED_EVAL_GRAD: 559 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 560 break; 561 // LCOV_EXCL_START 562 case CEED_EVAL_WEIGHT: 563 break; // Should not occur 564 case CEED_EVAL_DIV: 565 break; // TODO: Not implemented 566 case CEED_EVAL_CURL: 567 break; // TODO: Not implemented 568 // LCOV_EXCL_STOP 569 } 570 } 571 } else { 572 code << "\n // Note: Using full elements\n"; 573 code << " {\n"; 574 code << " // -- Input fields\n"; 575 for (CeedInt i = 0; i < num_input_fields; i++) { 576 code << " // ---- Input field " << i << "\n"; 577 code << " CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n"; 578 } 579 code << " // -- Output fields\n"; 580 for (CeedInt i = 0; i < num_output_fields; i++) { 581 code << " // ---- Output field " << i << "\n"; 582 code << " CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n"; 583 } 584 } 585 586 // Input and output buffers 587 code << "\n // -- QFunction inputs and outputs\n"; 588 code << " // ---- Inputs\n"; 589 code << " CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n"; 590 for (CeedInt i = 0; i < num_input_fields; i++) { 591 code << " // ------ Input field " << i << "\n"; 592 code << " inputs[" << i << "] = r_s_in_" << i << ";\n"; 593 } 594 code << " // ---- Outputs\n"; 595 code << " CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n"; 596 for (CeedInt i = 0; i < num_output_fields; i++) { 597 code << " // ------ Output field " << i << "\n"; 598 code << " outputs[" << i << "] = r_s_out_" << i << ";\n"; 599 } 600 601 // Apply QFunction 602 code << "\n // -- Apply QFunction\n"; 603 code << " " << qfunction_name << "(ctx, "; 604 if (dim != 3 || use_3d_slices) { 605 code << "1"; 606 } else { 607 code << "Q_1d"; 608 } 609 code << ", inputs, outputs);\n"; 610 611 // Copy or apply transpose grad, if needed 612 if (use_3d_slices) { 613 code << " // -- Output fields\n"; 614 for (CeedInt i = 0; i < num_output_fields; i++) { 615 std::string var_suffix = "_out_" + std::to_string(i); 616 std::string P_name = "P_1d" + var_suffix; 617 618 code << " // ---- Output field " << i << "\n"; 619 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 620 // Basis action 621 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 622 switch (eval_mode) { 623 case CEED_EVAL_NONE: 624 code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 625 code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 626 code << " }\n"; 627 break; // No action 628 case CEED_EVAL_INTERP: 629 code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 630 code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 631 code << " }\n"; 632 break; 633 case CEED_EVAL_GRAD: 634 code << " gradColloTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G" << var_suffix 635 << ", r_q" << var_suffix << ");\n"; 636 break; 637 // LCOV_EXCL_START 638 case CEED_EVAL_WEIGHT: 639 break; // Should not occur 640 case CEED_EVAL_DIV: 641 break; // TODO: Not implemented 642 case CEED_EVAL_CURL: 643 break; // TODO: Not implemented 644 // LCOV_EXCL_STOP 645 } 646 } 647 } 648 code << " }\n"; 649 return CEED_ERROR_SUCCESS; 650 } 651 652 //------------------------------------------------------------------------------ 653 // Build single operator kernel 654 //------------------------------------------------------------------------------ 655 extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op) { 656 bool is_tensor = true, use_3d_slices = false; 657 Ceed ceed; 658 CeedInt Q_1d, num_input_fields, num_output_fields, dim = 1; 659 CeedQFunctionField *qf_input_fields, *qf_output_fields; 660 CeedQFunction_Hip_gen *qf_data; 661 CeedQFunction qf; 662 CeedOperatorField *op_input_fields, *op_output_fields; 663 CeedOperator_Hip_gen *data; 664 std::ostringstream code; 665 666 { 667 bool is_setup_done; 668 669 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 670 if (is_setup_done) return CEED_ERROR_SUCCESS; 671 } 672 673 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 674 CeedCallBackend(CeedOperatorGetData(op, &data)); 675 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 676 CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 677 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 678 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 679 680 // Get operator data 681 CeedCallBackend(CeedOperatorBuildKernelData_Hip_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields, 682 qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices)); 683 if (dim == 0) dim = 1; 684 data->dim = dim; 685 if (Q_1d == 0) { 686 CeedInt Q; 687 688 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 689 Q_1d = Q; 690 } 691 data->Q_1d = Q_1d; 692 693 // Check for restriction only identity operator 694 { 695 bool is_identity_qf; 696 697 CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf)); 698 if (is_identity_qf) { 699 CeedEvalMode eval_mode_in, eval_mode_out; 700 701 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in)); 702 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out)); 703 CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND, 704 "Backend does not implement restriction only identity operators"); 705 } 706 } 707 708 // Load basis source files 709 // TODO: Add non-tensor, AtPoints 710 code << "// Tensor basis source\n"; 711 code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-templates.h>\n\n"; 712 code << "// CodeGen operator source\n"; 713 code << "#include <ceed/jit-source/hip/hip-gen-templates.h>\n\n"; 714 715 // Get QFunction name 716 std::string qfunction_name(qf_data->qfunction_name); 717 std::string operator_name; 718 719 operator_name = "CeedKernelHipGenOperator_" + qfunction_name; 720 721 // Define CEED_Q_VLA 722 code << "\n#undef CEED_Q_VLA\n"; 723 if (dim != 3 || use_3d_slices) { 724 code << "#define CEED_Q_VLA 1\n\n"; 725 } else { 726 code << "#define CEED_Q_VLA " << Q_1d << "\n\n"; 727 } 728 729 // Add user QFunction source 730 { 731 const char *source_path; 732 733 CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path)); 734 CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/hip/gen backend requires QFunction source code file"); 735 736 code << "// User QFunction source\n"; 737 code << "#include \"" << source_path << "\"\n\n"; 738 } 739 740 // Setup 741 code << "\n// -----------------------------------------------------------------------------\n"; 742 code << "// Operator Kernel\n"; 743 code << "// \n"; 744 code << "// d_[in,out]_i: CeedVector device array\n"; 745 code << "// r_[in,out]_e_i: Element vector register\n"; 746 code << "// r_[in,out]_q_i: Quadrature space vector register\n"; 747 code << "// r_[in,out]_s_i: Quadrature space slice vector register\n"; 748 code << "// \n"; 749 code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n"; 750 code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n"; 751 code << "// -----------------------------------------------------------------------------\n"; 752 code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n"; 753 code << "__global__ void " << operator_name 754 << "(CeedInt num_elem, void* ctx, FieldsInt_Hip indices, Fields_Hip fields, Fields_Hip B, Fields_Hip G, CeedScalar* W) {\n"; 755 756 // Scratch buffers 757 for (CeedInt i = 0; i < num_input_fields; i++) { 758 CeedEvalMode eval_mode; 759 760 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 761 if (eval_mode != CEED_EVAL_WEIGHT) { // Skip CEED_EVAL_WEIGHT 762 code << " const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n"; 763 } 764 } 765 for (CeedInt i = 0; i < num_output_fields; i++) { 766 code << " CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n"; 767 } 768 769 code << " const CeedInt dim = " << dim << ";\n"; 770 code << " const CeedInt Q_1d = " << Q_1d << ";\n"; 771 772 // Shared data 773 code << " extern __shared__ CeedScalar slice[];\n"; 774 code << " SharedData_Hip data;\n"; 775 code << " data.t_id_x = threadIdx.x;\n"; 776 code << " data.t_id_y = threadIdx.y;\n"; 777 code << " data.t_id_z = threadIdx.z;\n"; 778 code << " data.t_id = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n"; 779 code << " data.slice = slice + data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n"; 780 781 // Initialize constants, and matrices B and G 782 code << "\n // Input field constants and basis data\n"; 783 for (CeedInt i = 0; i < num_input_fields; i++) { 784 CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_input_fields[i], qf_input_fields[i], Q_1d, true, use_3d_slices)); 785 } 786 code << "\n // Output field constants and basis data\n"; 787 for (CeedInt i = 0; i < num_output_fields; i++) { 788 CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices)); 789 } 790 791 // Loop over all elements 792 code << "\n // Element loop\n"; 793 code << " __syncthreads();\n"; 794 code << " for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n"; 795 796 // -- Compute minimum buffer space needed 797 CeedInt max_rstr_buffer_size = 0; 798 799 for (CeedInt i = 0; i < num_input_fields; i++) { 800 CeedInt num_comp, elem_size; 801 CeedElemRestriction elem_rstr; 802 803 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 804 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 805 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 806 max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size); 807 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 808 } 809 for (CeedInt i = 0; i < num_output_fields; i++) { 810 CeedInt num_comp, elem_size; 811 CeedElemRestriction elem_rstr; 812 813 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 814 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 815 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 816 max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size); 817 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 818 } 819 code << " // Scratch restriction buffer space\n"; 820 code << " CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n"; 821 822 // -- Determine best input field processing order 823 CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX]; 824 825 for (CeedInt i = 0; i < num_input_fields; i++) { 826 field_rstr_in_buffer[i] = -1; 827 input_field_order[i] = -1; 828 } 829 { 830 bool is_ordered[CEED_FIELD_MAX]; 831 CeedInt curr_index = 0; 832 833 for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false; 834 for (CeedInt i = 0; i < num_input_fields; i++) { 835 CeedVector vec_i; 836 CeedElemRestriction rstr_i; 837 838 if (is_ordered[i]) continue; 839 field_rstr_in_buffer[i] = i; 840 is_ordered[i] = true; 841 input_field_order[curr_index] = i; 842 curr_index++; 843 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i)); 844 if (vec_i == CEED_VECTOR_NONE) continue; // CEED_EVAL_WEIGHT 845 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i)); 846 for (CeedInt j = i + 1; j < num_input_fields; j++) { 847 CeedVector vec_j; 848 CeedElemRestriction rstr_j; 849 850 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j)); 851 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j)); 852 if (rstr_i == rstr_j && vec_i == vec_j) { 853 field_rstr_in_buffer[j] = i; 854 is_ordered[j] = true; 855 input_field_order[curr_index] = j; 856 curr_index++; 857 } 858 } 859 } 860 } 861 862 // -- Input restriction and basis 863 code << " // -- Input field restrictions and basis actions\n"; 864 for (CeedInt i = 0; i < num_input_fields; i++) { 865 CeedInt f = input_field_order[i]; 866 867 code << " // ---- Input field " << f << "\n"; 868 869 // ---- Restriction 870 CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], Q_1d, 871 true, use_3d_slices)); 872 873 // ---- Basis action 874 CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, use_3d_slices)); 875 } 876 877 // -- Q function 878 CeedCallBackend(CeedOperatorBuildKernelQFunction_Hip_gen(code, data, dim, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, 879 op_output_fields, qf_output_fields, qfunction_name, Q_1d, use_3d_slices)); 880 881 // -- Output basis and restriction 882 code << "\n // -- Output field basis action and restrictions\n"; 883 for (CeedInt i = 0; i < num_output_fields; i++) { 884 code << " // ---- Output field " << i << "\n"; 885 886 // ---- Basis action 887 CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices)); 888 889 // ---- Restriction 890 CeedCallBackend( 891 CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices)); 892 } 893 894 // Close loop and function 895 code << " }\n"; 896 code << "}\n"; 897 code << "// -----------------------------------------------------------------------------\n\n"; 898 899 // View kernel for debugging 900 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated Operator Kernels:\n"); 901 CeedDebug(ceed, code.str().c_str()); 902 903 CeedInt block_sizes[3] = {0, 0, 0}; 904 CeedInt num_elem; 905 906 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 907 CeedCallBackend(BlockGridCalculate_Hip_gen(dim, num_elem, data->max_P_1d, Q_1d, block_sizes)); 908 CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE", 909 block_sizes[0] * block_sizes[1] * block_sizes[2])); 910 CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op)); 911 CeedCallBackend(CeedOperatorSetSetupDone(op)); 912 CeedCallBackend(CeedDestroy(&ceed)); 913 return CEED_ERROR_SUCCESS; 914 } 915 916 //------------------------------------------------------------------------------ 917