15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 37d8d0e25Snbeams // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 57d8d0e25Snbeams // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 73d576824SJeremy L Thompson 87d8d0e25Snbeams #define CEED_DEBUG_COLOR 12 97d8d0e25Snbeams 1049aac155SJeremy L Thompson #include <ceed.h> 11ec3da8bcSJed Brown #include <ceed/backend.h> 129e201c85SYohann #include <ceed/jit-tools.h> 132b730f8bSJeremy L Thompson 147d8d0e25Snbeams #include <iostream> 157d8d0e25Snbeams #include <sstream> 162b730f8bSJeremy L Thompson #include <string> 172b730f8bSJeremy L Thompson 180d0321e0SJeremy L Thompson #include "../hip-ref/ceed-hip-ref.h" 197d8d0e25Snbeams #include "../hip-shared/ceed-hip-shared.h" 20b2165e7aSSebastian Grimberg #include "../hip/ceed-hip-common.h" 217d8d0e25Snbeams #include "../hip/ceed-hip-compile.h" 222b730f8bSJeremy L Thompson #include "ceed-hip-gen.h" 23b3e1519bSnbeams 24b3e1519bSnbeams //------------------------------------------------------------------------------ 25b3e1519bSnbeams // Calculate the block size used for launching the operator kernel 26b3e1519bSnbeams //------------------------------------------------------------------------------ 272b730f8bSJeremy L Thompson extern "C" int BlockGridCalculate_Hip_gen(const CeedInt dim, const CeedInt num_elem, const CeedInt P_1d, const CeedInt Q_1d, CeedInt *block_sizes) { 289e201c85SYohann const CeedInt thread1d = CeedIntMax(Q_1d, P_1d); 29b3e1519bSnbeams if (dim == 1) { 309e201c85SYohann CeedInt elems_per_block = 64 * thread1d > 256 ? 256 / thread1d : 64; 31b7453713SJeremy L Thompson 329e201c85SYohann elems_per_block = elems_per_block > 0 ? elems_per_block : 1; 33b3e1519bSnbeams block_sizes[0] = thread1d; 34b3e1519bSnbeams block_sizes[1] = 1; 359e201c85SYohann block_sizes[2] = elems_per_block; 36b3e1519bSnbeams } else if (dim == 2) { 379e201c85SYohann const CeedInt elems_per_block = thread1d < 4 ? 16 : 2; 38b7453713SJeremy L Thompson 39b3e1519bSnbeams block_sizes[0] = thread1d; 40b3e1519bSnbeams block_sizes[1] = thread1d; 419e201c85SYohann block_sizes[2] = elems_per_block; 42b3e1519bSnbeams } else if (dim == 3) { 439e201c85SYohann const CeedInt elems_per_block = thread1d < 6 ? 4 : (thread1d < 8 ? 2 : 1); 44b7453713SJeremy L Thompson 45b3e1519bSnbeams block_sizes[0] = thread1d; 46b3e1519bSnbeams block_sizes[1] = thread1d; 479e201c85SYohann block_sizes[2] = elems_per_block; 48b3e1519bSnbeams } 49b3e1519bSnbeams return CEED_ERROR_SUCCESS; 50b3e1519bSnbeams } 51b3e1519bSnbeams 527d8d0e25Snbeams //------------------------------------------------------------------------------ 534b3e95d5SJeremy L Thompson // Determine type of operator 544b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 554b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Hip_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields, 564b3e95d5SJeremy L Thompson CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields, 574b3e95d5SJeremy L Thompson CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor, 584b3e95d5SJeremy L Thompson bool *use_3d_slices) { 594b3e95d5SJeremy L Thompson // Find dim, P_1d, Q_1d 604b3e95d5SJeremy L Thompson *max_P_1d = 0; 614b3e95d5SJeremy L Thompson *Q_1d = 0; 624b3e95d5SJeremy L Thompson *dim = 0; 634b3e95d5SJeremy L Thompson *is_tensor = true; 644b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 654b3e95d5SJeremy L Thompson CeedBasis basis; 664b3e95d5SJeremy L Thompson 674b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 684b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 694b3e95d5SJeremy L Thompson bool is_field_tensor; 704b3e95d5SJeremy L Thompson CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0; 714b3e95d5SJeremy L Thompson 724b3e95d5SJeremy L Thompson // Collect dim, P_1d, and Q_1d 734b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 744b3e95d5SJeremy L Thompson CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 754b3e95d5SJeremy L Thompson *is_tensor = *is_tensor && is_field_tensor; 764b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 774b3e95d5SJeremy L Thompson *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 784b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 794b3e95d5SJeremy L Thompson CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 804b3e95d5SJeremy L Thompson *dim = field_dim; 814b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 824b3e95d5SJeremy L Thompson CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 834b3e95d5SJeremy L Thompson *Q_1d = field_Q_1d; 844b3e95d5SJeremy L Thompson } 854b3e95d5SJeremy L Thompson } 864b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 874b3e95d5SJeremy L Thompson CeedBasis basis; 884b3e95d5SJeremy L Thompson 894b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 904b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 914b3e95d5SJeremy L Thompson bool is_field_tensor; 924b3e95d5SJeremy L Thompson CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0; 934b3e95d5SJeremy L Thompson 944b3e95d5SJeremy L Thompson // Collect dim, P_1d, and Q_1d 954b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 964b3e95d5SJeremy L Thompson CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 974b3e95d5SJeremy L Thompson *is_tensor = *is_tensor && is_field_tensor; 984b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 994b3e95d5SJeremy L Thompson *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 1004b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 1014b3e95d5SJeremy L Thompson CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 1024b3e95d5SJeremy L Thompson *dim = field_dim; 1034b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 1044b3e95d5SJeremy L Thompson CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 1054b3e95d5SJeremy L Thompson *Q_1d = field_Q_1d; 1064b3e95d5SJeremy L Thompson } 1074b3e95d5SJeremy L Thompson } 1084b3e95d5SJeremy L Thompson 1094b3e95d5SJeremy L Thompson // Only use 3D collocated gradient parallelization strategy when gradient is computed 1104b3e95d5SJeremy L Thompson *use_3d_slices = false; 1114b3e95d5SJeremy L Thompson if (*dim == 3) { 1124b3e95d5SJeremy L Thompson bool was_grad_found = false; 1134b3e95d5SJeremy L Thompson 1144b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 1154b3e95d5SJeremy L Thompson CeedEvalMode eval_mode; 1164b3e95d5SJeremy L Thompson 1174b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 1184b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_GRAD) { 1194b3e95d5SJeremy L Thompson CeedBasis_Hip_shared *basis_data; 1204b3e95d5SJeremy L Thompson CeedBasis basis; 1214b3e95d5SJeremy L Thompson 1224b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 1234b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 1244b3e95d5SJeremy L Thompson *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 1254b3e95d5SJeremy L Thompson was_grad_found = true; 1264b3e95d5SJeremy L Thompson } 1274b3e95d5SJeremy L Thompson } 1284b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 1294b3e95d5SJeremy L Thompson CeedEvalMode eval_mode; 1304b3e95d5SJeremy L Thompson 1314b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1324b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_GRAD) { 1334b3e95d5SJeremy L Thompson CeedBasis_Hip_shared *basis_data; 1344b3e95d5SJeremy L Thompson CeedBasis basis; 1354b3e95d5SJeremy L Thompson 1364b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 1374b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 1384b3e95d5SJeremy L Thompson *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 1394b3e95d5SJeremy L Thompson was_grad_found = true; 1404b3e95d5SJeremy L Thompson } 1414b3e95d5SJeremy L Thompson } 1424b3e95d5SJeremy L Thompson } 1434b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 1444b3e95d5SJeremy L Thompson } 1454b3e95d5SJeremy L Thompson 1464b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 1474b3e95d5SJeremy L Thompson // Setup fields 1484b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 1494b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedOperatorField op_field, 1504b3e95d5SJeremy L Thompson CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool use_3d_slices) { 1514b3e95d5SJeremy L Thompson std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 1524b3e95d5SJeremy L Thompson std::string P_name = "P_1d" + var_suffix, Q_name = "Q_1d"; 1534b3e95d5SJeremy L Thompson std::string option_name = (is_input ? "inputs" : "outputs"); 1544b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 1554b3e95d5SJeremy L Thompson CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 1564b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 1574b3e95d5SJeremy L Thompson CeedBasis_Hip_shared *basis_data; 1584b3e95d5SJeremy L Thompson CeedBasis basis; 1594b3e95d5SJeremy L Thompson 1604b3e95d5SJeremy L Thompson code << " // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n"; 1614b3e95d5SJeremy L Thompson 1624b3e95d5SJeremy L Thompson // Get field data 1634b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 1644b3e95d5SJeremy L Thompson if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 1654b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1664b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 1674b3e95d5SJeremy L Thompson } 1684b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 1694b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 1704b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 1714b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 1724b3e95d5SJeremy L Thompson } 1734b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 1744b3e95d5SJeremy L Thompson 1754b3e95d5SJeremy L Thompson // Set field constants 1764b3e95d5SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) { 1774b3e95d5SJeremy L Thompson code << " const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n"; 1784b3e95d5SJeremy L Thompson code << " const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n"; 1794b3e95d5SJeremy L Thompson } 1804b3e95d5SJeremy L Thompson 1814b3e95d5SJeremy L Thompson // Load basis data 1824b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 1834b3e95d5SJeremy L Thompson switch (eval_mode) { 1844b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 1854b3e95d5SJeremy L Thompson break; 1864b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 1874b3e95d5SJeremy L Thompson if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 1884b3e95d5SJeremy L Thompson else data->B.outputs[i] = basis_data->d_interp_1d; 1894b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n"; 1904b3e95d5SJeremy L Thompson code << " loadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 1914b3e95d5SJeremy L Thompson break; 1924b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 1934b3e95d5SJeremy L Thompson if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 1944b3e95d5SJeremy L Thompson else data->B.outputs[i] = basis_data->d_interp_1d; 1954b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n"; 1964b3e95d5SJeremy L Thompson code << " loadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 1974b3e95d5SJeremy L Thompson if (use_3d_slices) { 1984b3e95d5SJeremy L Thompson if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d; 1994b3e95d5SJeremy L Thompson else data->G.outputs[i] = basis_data->d_collo_grad_1d; 2004b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n"; 2014b3e95d5SJeremy L Thompson code << " loadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 2024b3e95d5SJeremy L Thompson } else { 2034b3e95d5SJeremy L Thompson bool has_collo_grad = basis_data->d_collo_grad_1d; 2044b3e95d5SJeremy L Thompson 2054b3e95d5SJeremy L Thompson if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 2064b3e95d5SJeremy L Thompson else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 2074b3e95d5SJeremy L Thompson if (has_collo_grad) { 2084b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n"; 2094b3e95d5SJeremy L Thompson code << " loadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 2104b3e95d5SJeremy L Thompson } else { 2114b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * P_1d << "];\n"; 2124b3e95d5SJeremy L Thompson code << " loadMatrix<" << P_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 2134b3e95d5SJeremy L Thompson } 2144b3e95d5SJeremy L Thompson } 2154b3e95d5SJeremy L Thompson break; 2164b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 2174b3e95d5SJeremy L Thompson break; // No action 2184b3e95d5SJeremy L Thompson // LCOV_EXCL_START 2194b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 2204b3e95d5SJeremy L Thompson break; // TODO: Not implemented 2214b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 2224b3e95d5SJeremy L Thompson break; // TODO: Not implemented 2234b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 2244b3e95d5SJeremy L Thompson } 2254b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 2264b3e95d5SJeremy L Thompson } 2274b3e95d5SJeremy L Thompson 2284b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 2294b3e95d5SJeremy L Thompson // Restriction 2304b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 2314b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim, 232e93651e5SJeremy L Thompson CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field, 233e93651e5SJeremy L Thompson CeedInt Q_1d, bool is_input, bool use_3d_slices) { 2344b3e95d5SJeremy L Thompson std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 2354b3e95d5SJeremy L Thompson std::string P_name = "P_1d" + var_suffix; 2364b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 2374b3e95d5SJeremy L Thompson CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 2384b3e95d5SJeremy L Thompson CeedSize l_size; 2394b3e95d5SJeremy L Thompson CeedElemRestriction_Hip *rstr_data; 2404b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 2414b3e95d5SJeremy L Thompson CeedBasis basis; 2424b3e95d5SJeremy L Thompson 2434b3e95d5SJeremy L Thompson // Get field data 2444b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 2454b3e95d5SJeremy L Thompson if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 2464b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 2474b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 2484b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 2494b3e95d5SJeremy L Thompson } 2504b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 2514b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 2524b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 2534b3e95d5SJeremy L Thompson } 2544b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 2554b3e95d5SJeremy L Thompson 2564b3e95d5SJeremy L Thompson // Restriction 2574b3e95d5SJeremy L Thompson if (is_input) { 2584b3e95d5SJeremy L Thompson // Input 259e93651e5SJeremy L Thompson // Input 260e93651e5SJeremy L Thompson if (field_input_buffer[i] != i) { 261e93651e5SJeremy L Thompson std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]); 262e93651e5SJeremy L Thompson 263e93651e5SJeremy L Thompson // Restriction was already done for previous input 264e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n"; 265e93651e5SJeremy L Thompson } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices)) { 2664b3e95d5SJeremy L Thompson bool is_strided; 2674b3e95d5SJeremy L Thompson 268e93651e5SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) { 269e93651e5SJeremy L Thompson // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated 2704b3e95d5SJeremy L Thompson code << " CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n"; 271e93651e5SJeremy L Thompson } else { 272e93651e5SJeremy L Thompson // Otherwise we're using the scratch space 273e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 274e93651e5SJeremy L Thompson } 2754b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 2764b3e95d5SJeremy L Thompson if (!is_strided) { 2774b3e95d5SJeremy L Thompson CeedInt comp_stride; 2784b3e95d5SJeremy L Thompson 2794b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 2804b3e95d5SJeremy L Thompson code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 2814b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 2824b3e95d5SJeremy L Thompson code << " // CompStride: " << comp_stride << "\n"; 2834b3e95d5SJeremy L Thompson data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 2844b3e95d5SJeremy L Thompson code << " readDofsOffset" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" << var_suffix 2854b3e95d5SJeremy L Thompson << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n"; 2864b3e95d5SJeremy L Thompson } else { 2874b3e95d5SJeremy L Thompson bool has_backend_strides; 2884b3e95d5SJeremy L Thompson CeedInt num_elem; 2894b3e95d5SJeremy L Thompson 2904b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 2914b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 2924b3e95d5SJeremy L Thompson CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 2934b3e95d5SJeremy L Thompson 2944b3e95d5SJeremy L Thompson if (!has_backend_strides) { 2954b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 2964b3e95d5SJeremy L Thompson } 2974b3e95d5SJeremy L Thompson code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 2984b3e95d5SJeremy L Thompson code << " readDofsStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << "," 2994b3e95d5SJeremy L Thompson << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n"; 3004b3e95d5SJeremy L Thompson } 3014b3e95d5SJeremy L Thompson } 3024b3e95d5SJeremy L Thompson } else { 3034b3e95d5SJeremy L Thompson // Output 3044b3e95d5SJeremy L Thompson bool is_strided; 3054b3e95d5SJeremy L Thompson 3064b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 3074b3e95d5SJeremy L Thompson if (!is_strided) { 3084b3e95d5SJeremy L Thompson CeedInt comp_stride; 3094b3e95d5SJeremy L Thompson 3104b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 3114b3e95d5SJeremy L Thompson code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 3124b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 3134b3e95d5SJeremy L Thompson code << " // CompStride: " << comp_stride << "\n"; 3144b3e95d5SJeremy L Thompson data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets; 3154b3e95d5SJeremy L Thompson code << " writeDofsOffset" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" << var_suffix 3164b3e95d5SJeremy L Thompson << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n"; 3174b3e95d5SJeremy L Thompson } else { 3184b3e95d5SJeremy L Thompson bool has_backend_strides; 3194b3e95d5SJeremy L Thompson CeedInt num_elem; 3204b3e95d5SJeremy L Thompson 3214b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 3224b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 3234b3e95d5SJeremy L Thompson CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 3244b3e95d5SJeremy L Thompson 3254b3e95d5SJeremy L Thompson if (!has_backend_strides) { 3264b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 3274b3e95d5SJeremy L Thompson } 3284b3e95d5SJeremy L Thompson code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 3294b3e95d5SJeremy L Thompson code << " writeDofsStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << "," 3304b3e95d5SJeremy L Thompson << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n"; 3314b3e95d5SJeremy L Thompson } 3324b3e95d5SJeremy L Thompson } 3334b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 3344b3e95d5SJeremy L Thompson } 3354b3e95d5SJeremy L Thompson 3364b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 3374b3e95d5SJeremy L Thompson // Basis 3384b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 3394b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt dim, 3404b3e95d5SJeremy L Thompson CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, 3414b3e95d5SJeremy L Thompson bool use_3d_slices) { 3424b3e95d5SJeremy L Thompson std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 3434b3e95d5SJeremy L Thompson std::string P_name = "P_1d" + var_suffix, Q_name = "Q_1d"; 3444b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 3454b3e95d5SJeremy L Thompson CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 3464b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 3474b3e95d5SJeremy L Thompson CeedBasis basis; 3484b3e95d5SJeremy L Thompson 3494b3e95d5SJeremy L Thompson // Get field data 3504b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 3514b3e95d5SJeremy L Thompson if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 3524b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 3534b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 3544b3e95d5SJeremy L Thompson } 3554b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 3564b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 3574b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 3584b3e95d5SJeremy L Thompson } 3594b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 3604b3e95d5SJeremy L Thompson 3614b3e95d5SJeremy L Thompson // Basis 3624b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 3634b3e95d5SJeremy L Thompson if (is_input) { 3644b3e95d5SJeremy L Thompson switch (eval_mode) { 3654b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 3664b3e95d5SJeremy L Thompson if (!use_3d_slices) { 3674b3e95d5SJeremy L Thompson code << " CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n"; 3684b3e95d5SJeremy L Thompson } 3694b3e95d5SJeremy L Thompson break; 3704b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 3714b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 3724b3e95d5SJeremy L Thompson code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name 3734b3e95d5SJeremy L Thompson << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n"; 3744b3e95d5SJeremy L Thompson break; 3754b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 3764b3e95d5SJeremy L Thompson if (use_3d_slices) { 3774b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 3784b3e95d5SJeremy L Thompson code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name 3794b3e95d5SJeremy L Thompson << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n"; 3804b3e95d5SJeremy L Thompson } else { 3814b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n"; 3824b3e95d5SJeremy L Thompson code << " Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" << var_suffix 3834b3e95d5SJeremy L Thompson << ", P_1d" << var_suffix << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q" 3844b3e95d5SJeremy L Thompson << var_suffix << ");\n"; 3854b3e95d5SJeremy L Thompson } 3864b3e95d5SJeremy L Thompson break; 3874b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: { 3884b3e95d5SJeremy L Thompson CeedBasis_Hip_shared *basis_data; 3894b3e95d5SJeremy L Thompson 3904b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[" << Q_name << "];\n"; 3914b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 3924b3e95d5SJeremy L Thompson data->W = basis_data->d_q_weight_1d; 3934b3e95d5SJeremy L Thompson code << " Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n"; 3944b3e95d5SJeremy L Thompson break; 3954b3e95d5SJeremy L Thompson } 3964b3e95d5SJeremy L Thompson // LCOV_EXCL_START 3974b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 3984b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 3994b3e95d5SJeremy L Thompson break; // TODO: Not implemented 4004b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 4014b3e95d5SJeremy L Thompson } 4024b3e95d5SJeremy L Thompson } else { 4034b3e95d5SJeremy L Thompson switch (eval_mode) { 4044b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 4054b3e95d5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n"; 4064b3e95d5SJeremy L Thompson break; // No action 4074b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 408e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 4094b3e95d5SJeremy L Thompson code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name 4104b3e95d5SJeremy L Thompson << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 4114b3e95d5SJeremy L Thompson break; 4124b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 413e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 4144b3e95d5SJeremy L Thompson if (use_3d_slices) { 4154b3e95d5SJeremy L Thompson code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name 4164b3e95d5SJeremy L Thompson << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 4174b3e95d5SJeremy L Thompson } else { 4184b3e95d5SJeremy L Thompson code << " GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" 4194b3e95d5SJeremy L Thompson << var_suffix << ", " << P_name << "," << Q_name << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix 4204b3e95d5SJeremy L Thompson << ", r_e" << var_suffix << ");\n"; 4214b3e95d5SJeremy L Thompson } 4224b3e95d5SJeremy L Thompson break; 4234b3e95d5SJeremy L Thompson // LCOV_EXCL_START 4244b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 4254b3e95d5SJeremy L Thompson break; // Should not occur 4264b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 4274b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 4284b3e95d5SJeremy L Thompson break; // TODO: Not implemented 4294b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 4304b3e95d5SJeremy L Thompson } 4314b3e95d5SJeremy L Thompson } 4324b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 4334b3e95d5SJeremy L Thompson } 4344b3e95d5SJeremy L Thompson 4354b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 4364b3e95d5SJeremy L Thompson // QFunction 4374b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 4384b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt dim, CeedInt num_input_fields, 4394b3e95d5SJeremy L Thompson CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields, 4404b3e95d5SJeremy L Thompson CeedInt num_output_fields, CeedOperatorField *op_output_fields, 4414b3e95d5SJeremy L Thompson CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d, 4424b3e95d5SJeremy L Thompson bool use_3d_slices) { 4434b3e95d5SJeremy L Thompson std::string Q_name = "Q_1d"; 4444b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 4454b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 4464b3e95d5SJeremy L Thompson 4474b3e95d5SJeremy L Thompson // Setup output arays 4484b3e95d5SJeremy L Thompson code << "\n // -- Output field setup\n"; 4494b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 4504b3e95d5SJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 4514b3e95d5SJeremy L Thompson 4524b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 4534b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 4544b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) { 4554b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 4564b3e95d5SJeremy L Thompson } 4574b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_GRAD) { 4584b3e95d5SJeremy L Thompson if (use_3d_slices) { 4594b3e95d5SJeremy L Thompson // Accumulator for gradient slices 4604b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 4614b3e95d5SJeremy L Thompson code << " for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n"; 4624b3e95d5SJeremy L Thompson code << " r_q" << var_suffix << "[i] = 0.0;\n"; 4634b3e95d5SJeremy L Thompson code << " }\n"; 4644b3e95d5SJeremy L Thompson } else { 4654b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n"; 4664b3e95d5SJeremy L Thompson } 4674b3e95d5SJeremy L Thompson } 4684b3e95d5SJeremy L Thompson } 4694b3e95d5SJeremy L Thompson 4704b3e95d5SJeremy L Thompson // We treat quadrature points per slice in 3d to save registers 4714b3e95d5SJeremy L Thompson if (use_3d_slices) { 4724b3e95d5SJeremy L Thompson code << "\n // Note: Using planes of 3D elements\n"; 4734b3e95d5SJeremy L Thompson code << " #pragma unroll\n"; 4744b3e95d5SJeremy L Thompson code << " for (CeedInt q = 0; q < " << Q_name << "; q++) {\n"; 4754b3e95d5SJeremy L Thompson code << " // -- Input fields\n"; 4764b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 4774b3e95d5SJeremy L Thompson std::string var_suffix = "_in_" + std::to_string(i); 4784b3e95d5SJeremy L Thompson 4794b3e95d5SJeremy L Thompson code << " // ---- Input field " << i << "\n"; 4804b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 4814b3e95d5SJeremy L Thompson // Basis action 4824b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 4834b3e95d5SJeremy L Thompson switch (eval_mode) { 4844b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 4854b3e95d5SJeremy L Thompson bool is_strided; 4864b3e95d5SJeremy L Thompson 4874b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 4884b3e95d5SJeremy L Thompson 4894b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 4904b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 4914b3e95d5SJeremy L Thompson if (is_strided) { 4924b3e95d5SJeremy L Thompson bool has_backend_strides; 4934b3e95d5SJeremy L Thompson CeedInt num_elem, elem_size; 4944b3e95d5SJeremy L Thompson 4954b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 4964b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 4974b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 4984b3e95d5SJeremy L Thompson CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 4994b3e95d5SJeremy L Thompson 5004b3e95d5SJeremy L Thompson if (!has_backend_strides) { 5014b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 5024b3e95d5SJeremy L Thompson } 5034b3e95d5SJeremy L Thompson code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 5044b3e95d5SJeremy L Thompson code << " readSliceQuadsStrided3d<num_comp" << var_suffix << ", " << Q_name << "," << strides[0] << "," << strides[1] << "," 5054b3e95d5SJeremy L Thompson << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n"; 5064b3e95d5SJeremy L Thompson } else { 5074b3e95d5SJeremy L Thompson CeedSize l_size = 0; 5084b3e95d5SJeremy L Thompson CeedInt comp_stride; 5094b3e95d5SJeremy L Thompson CeedElemRestriction_Hip *rstr_data; 5104b3e95d5SJeremy L Thompson 5114b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 5124b3e95d5SJeremy L Thompson code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 5134b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 5144b3e95d5SJeremy L Thompson code << " // CompStride: " << comp_stride << "\n"; 5154b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 5164b3e95d5SJeremy L Thompson data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 5174b3e95d5SJeremy L Thompson code << " readSliceQuadsOffset3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix 5184b3e95d5SJeremy L Thompson << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n"; 5194b3e95d5SJeremy L Thompson } 5204b3e95d5SJeremy L Thompson break; 5214b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 5224b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 5234b3e95d5SJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n"; 5244b3e95d5SJeremy L Thompson code << " r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n"; 5254b3e95d5SJeremy L Thompson code << " }\n"; 5264b3e95d5SJeremy L Thompson break; 5274b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 5284b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 5294b3e95d5SJeremy L Thompson code << " gradCollo3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix << ", r_s" 5304b3e95d5SJeremy L Thompson << var_suffix << ");\n"; 5314b3e95d5SJeremy L Thompson break; 5324b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 5334b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[1];\n"; 5344b3e95d5SJeremy L Thompson code << " r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n"; 5354b3e95d5SJeremy L Thompson break; // No action 5364b3e95d5SJeremy L Thompson // LCOV_EXCL_START 5374b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 5384b3e95d5SJeremy L Thompson break; // TODO: Not implemented 5394b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 5404b3e95d5SJeremy L Thompson break; // TODO: Not implemented 5414b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 5424b3e95d5SJeremy L Thompson } 5434b3e95d5SJeremy L Thompson } 5444b3e95d5SJeremy L Thompson code << "\n // -- Output fields\n"; 5454b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 5464b3e95d5SJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 5474b3e95d5SJeremy L Thompson 5484b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 5494b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 5504b3e95d5SJeremy L Thompson // Basis action 5514b3e95d5SJeremy L Thompson switch (eval_mode) { 5524b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 5534b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 5544b3e95d5SJeremy L Thompson break; // No action 5554b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 5564b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 5574b3e95d5SJeremy L Thompson break; 5584b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 5594b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 5604b3e95d5SJeremy L Thompson break; 5614b3e95d5SJeremy L Thompson // LCOV_EXCL_START 5624b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 5634b3e95d5SJeremy L Thompson break; // Should not occur 5644b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 5654b3e95d5SJeremy L Thompson break; // TODO: Not implemented 5664b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 5674b3e95d5SJeremy L Thompson break; // TODO: Not implemented 5684b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 5694b3e95d5SJeremy L Thompson } 5704b3e95d5SJeremy L Thompson } 5714b3e95d5SJeremy L Thompson } else { 5724b3e95d5SJeremy L Thompson code << "\n // Note: Using full elements\n"; 5734b3e95d5SJeremy L Thompson code << " {\n"; 5744b3e95d5SJeremy L Thompson code << " // -- Input fields\n"; 5754b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 5764b3e95d5SJeremy L Thompson code << " // ---- Input field " << i << "\n"; 5774b3e95d5SJeremy L Thompson code << " CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n"; 5784b3e95d5SJeremy L Thompson } 5794b3e95d5SJeremy L Thompson code << " // -- Output fields\n"; 5804b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 5814b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 5824b3e95d5SJeremy L Thompson code << " CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n"; 5834b3e95d5SJeremy L Thompson } 5844b3e95d5SJeremy L Thompson } 5854b3e95d5SJeremy L Thompson 5864b3e95d5SJeremy L Thompson // Input and output buffers 5874b3e95d5SJeremy L Thompson code << "\n // -- QFunction inputs and outputs\n"; 5884b3e95d5SJeremy L Thompson code << " // ---- Inputs\n"; 5894b3e95d5SJeremy L Thompson code << " CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n"; 5904b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 5914b3e95d5SJeremy L Thompson code << " // ------ Input field " << i << "\n"; 5924b3e95d5SJeremy L Thompson code << " inputs[" << i << "] = r_s_in_" << i << ";\n"; 5934b3e95d5SJeremy L Thompson } 5944b3e95d5SJeremy L Thompson code << " // ---- Outputs\n"; 5954b3e95d5SJeremy L Thompson code << " CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n"; 5964b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 5974b3e95d5SJeremy L Thompson code << " // ------ Output field " << i << "\n"; 5984b3e95d5SJeremy L Thompson code << " outputs[" << i << "] = r_s_out_" << i << ";\n"; 5994b3e95d5SJeremy L Thompson } 6004b3e95d5SJeremy L Thompson 6014b3e95d5SJeremy L Thompson // Apply QFunction 6024b3e95d5SJeremy L Thompson code << "\n // -- Apply QFunction\n"; 6034b3e95d5SJeremy L Thompson code << " " << qfunction_name << "(ctx, "; 6044b3e95d5SJeremy L Thompson if (dim != 3 || use_3d_slices) { 6054b3e95d5SJeremy L Thompson code << "1"; 6064b3e95d5SJeremy L Thompson } else { 6074b3e95d5SJeremy L Thompson code << "Q_1d"; 6084b3e95d5SJeremy L Thompson } 6094b3e95d5SJeremy L Thompson code << ", inputs, outputs);\n"; 6104b3e95d5SJeremy L Thompson 6114b3e95d5SJeremy L Thompson // Copy or apply transpose grad, if needed 6124b3e95d5SJeremy L Thompson if (use_3d_slices) { 6134b3e95d5SJeremy L Thompson code << " // -- Output fields\n"; 6144b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 6154b3e95d5SJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 6164b3e95d5SJeremy L Thompson std::string P_name = "P_1d" + var_suffix; 6174b3e95d5SJeremy L Thompson 6184b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 6194b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 6204b3e95d5SJeremy L Thompson // Basis action 6214b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 6224b3e95d5SJeremy L Thompson switch (eval_mode) { 6234b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 6244b3e95d5SJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 6254b3e95d5SJeremy L Thompson code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 6264b3e95d5SJeremy L Thompson code << " }\n"; 6274b3e95d5SJeremy L Thompson break; // No action 6284b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 6294b3e95d5SJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 6304b3e95d5SJeremy L Thompson code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 6314b3e95d5SJeremy L Thompson code << " }\n"; 6324b3e95d5SJeremy L Thompson break; 6334b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 6344b3e95d5SJeremy L Thompson code << " gradColloTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G" << var_suffix 6354b3e95d5SJeremy L Thompson << ", r_q" << var_suffix << ");\n"; 6364b3e95d5SJeremy L Thompson break; 6374b3e95d5SJeremy L Thompson // LCOV_EXCL_START 6384b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 6394b3e95d5SJeremy L Thompson break; // Should not occur 6404b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 6414b3e95d5SJeremy L Thompson break; // TODO: Not implemented 6424b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 6434b3e95d5SJeremy L Thompson break; // TODO: Not implemented 6444b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 6454b3e95d5SJeremy L Thompson } 6464b3e95d5SJeremy L Thompson } 6474b3e95d5SJeremy L Thompson } 6484b3e95d5SJeremy L Thompson code << " }\n"; 6494b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 6504b3e95d5SJeremy L Thompson } 6514b3e95d5SJeremy L Thompson 6524b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 6539e201c85SYohann // Build single operator kernel 6547d8d0e25Snbeams //------------------------------------------------------------------------------ 655eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op) { 6564b3e95d5SJeremy L Thompson bool is_tensor = true, use_3d_slices = false; 6577d8d0e25Snbeams Ceed ceed; 6584b3e95d5SJeremy L Thompson CeedInt Q_1d, num_input_fields, num_output_fields, dim = 1; 659b7453713SJeremy L Thompson CeedQFunctionField *qf_input_fields, *qf_output_fields; 660b7453713SJeremy L Thompson CeedQFunction_Hip_gen *qf_data; 661b7453713SJeremy L Thompson CeedQFunction qf; 662b7453713SJeremy L Thompson CeedOperatorField *op_input_fields, *op_output_fields; 663b7453713SJeremy L Thompson CeedOperator_Hip_gen *data; 6644b3e95d5SJeremy L Thompson std::ostringstream code; 6654b3e95d5SJeremy L Thompson 6664b3e95d5SJeremy L Thompson { 6674b3e95d5SJeremy L Thompson bool is_setup_done; 668b7453713SJeremy L Thompson 669b7453713SJeremy L Thompson CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 670b7453713SJeremy L Thompson if (is_setup_done) return CEED_ERROR_SUCCESS; 6714b3e95d5SJeremy L Thompson } 672b7453713SJeremy L Thompson 673b7453713SJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 674b7453713SJeremy L Thompson CeedCallBackend(CeedOperatorGetData(op, &data)); 675b7453713SJeremy L Thompson CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 676b7453713SJeremy L Thompson CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 677b7453713SJeremy L Thompson CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 678b7453713SJeremy L Thompson CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 6797d8d0e25Snbeams 6804b3e95d5SJeremy L Thompson // Get operator data 6814b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelData_Hip_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields, 6824b3e95d5SJeremy L Thompson qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices)); 6834b3e95d5SJeremy L Thompson if (dim == 0) dim = 1; 6844b3e95d5SJeremy L Thompson data->dim = dim; 6854b3e95d5SJeremy L Thompson if (Q_1d == 0) { 6864b3e95d5SJeremy L Thompson CeedInt Q; 6874b3e95d5SJeremy L Thompson 6884b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 6894b3e95d5SJeremy L Thompson Q_1d = Q; 6904b3e95d5SJeremy L Thompson } 6914b3e95d5SJeremy L Thompson data->Q_1d = Q_1d; 6924b3e95d5SJeremy L Thompson 6930b454692Sjeremylt // Check for restriction only identity operator 6944b3e95d5SJeremy L Thompson { 6954b3e95d5SJeremy L Thompson bool is_identity_qf; 6964b3e95d5SJeremy L Thompson 6972b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf)); 6980b454692Sjeremylt if (is_identity_qf) { 6999e201c85SYohann CeedEvalMode eval_mode_in, eval_mode_out; 700b7453713SJeremy L Thompson 7012b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in)); 7022b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out)); 7036574a04fSJeremy L Thompson CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND, 7046574a04fSJeremy L Thompson "Backend does not implement restriction only identity operators"); 7050b454692Sjeremylt } 7064b3e95d5SJeremy L Thompson } 707b2165e7aSSebastian Grimberg 708b2165e7aSSebastian Grimberg // Load basis source files 7094b3e95d5SJeremy L Thompson // TODO: Add non-tensor, AtPoints 7109e201c85SYohann { 71122070f95SJeremy L Thompson char *tensor_basis_kernel_source; 71222070f95SJeremy L Thompson const char *tensor_basis_kernel_path; 713b7453713SJeremy L Thompson 7142b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-shared-basis-tensor-templates.h", &tensor_basis_kernel_path)); 71523d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Tensor Basis Kernel Source -----\n"); 7162b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source)); 7179e201c85SYohann code << tensor_basis_kernel_source; 7182b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&tensor_basis_kernel_path)); 7192b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&tensor_basis_kernel_source)); 7209e201c85SYohann } 7219e201c85SYohann { 72222070f95SJeremy L Thompson char *hip_gen_template_source; 72322070f95SJeremy L Thompson const char *hip_gen_template_path; 724b7453713SJeremy L Thompson 7252b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-gen-templates.h", &hip_gen_template_path)); 72623d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Hip-Gen Template Source -----\n"); 7272b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, hip_gen_template_path, &hip_gen_template_source)); 7289e201c85SYohann code << hip_gen_template_source; 7292b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&hip_gen_template_path)); 7302b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&hip_gen_template_source)); 7319e201c85SYohann } 7327d8d0e25Snbeams 7334b3e95d5SJeremy L Thompson // Get QFunction name 7344b3e95d5SJeremy L Thompson std::string qfunction_name(qf_data->qfunction_name); 7354b3e95d5SJeremy L Thompson std::string operator_name; 7364b3e95d5SJeremy L Thompson 73709095acaSJeremy L Thompson operator_name = "CeedKernelHipGenOperator_" + qfunction_name; 7387d8d0e25Snbeams 7399e201c85SYohann // Define CEED_Q_VLA 7409e201c85SYohann code << "\n#undef CEED_Q_VLA\n"; 7414b3e95d5SJeremy L Thompson if (dim != 3 || use_3d_slices) { 7429e201c85SYohann code << "#define CEED_Q_VLA 1\n\n"; 7439e201c85SYohann } else { 7449e201c85SYohann code << "#define CEED_Q_VLA " << Q_1d << "\n\n"; 7459e201c85SYohann } 7469e201c85SYohann 7474b3e95d5SJeremy L Thompson // Add user QFunction source 7484b3e95d5SJeremy L Thompson { 7494b3e95d5SJeremy L Thompson std::string qfunction_source(qf_data->qfunction_source); 7504b3e95d5SJeremy L Thompson 75109095acaSJeremy L Thompson code << qfunction_source; 7524b3e95d5SJeremy L Thompson } 7537d8d0e25Snbeams 7547d8d0e25Snbeams // Setup 7557d8d0e25Snbeams code << "\n// -----------------------------------------------------------------------------\n"; 7564b3e95d5SJeremy L Thompson code << "// Operator Kernel\n"; 7574b3e95d5SJeremy L Thompson code << "// \n"; 7584b3e95d5SJeremy L Thompson code << "// d_[in,out]_i: CeedVector device array\n"; 7594b3e95d5SJeremy L Thompson code << "// r_[in,out]_e_i: Element vector register\n"; 7604b3e95d5SJeremy L Thompson code << "// r_[in,out]_q_i: Quadrature space vector register\n"; 7614b3e95d5SJeremy L Thompson code << "// r_[in,out]_s_i: Quadrature space slice vector register\n"; 7624b3e95d5SJeremy L Thompson code << "// \n"; 7634b3e95d5SJeremy L Thompson code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n"; 7644b3e95d5SJeremy L Thompson code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n"; 7654b3e95d5SJeremy L Thompson code << "// -----------------------------------------------------------------------------\n"; 766b3e1519bSnbeams code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n"; 7672b730f8bSJeremy L Thompson code << "__global__ void " << operator_name 7682b730f8bSJeremy L Thompson << "(CeedInt num_elem, void* ctx, FieldsInt_Hip indices, Fields_Hip fields, Fields_Hip B, Fields_Hip G, CeedScalar* W) {\n"; 7694b3e95d5SJeremy L Thompson 7704b3e95d5SJeremy L Thompson // Scratch buffers 7719e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 7724b3e95d5SJeremy L Thompson CeedEvalMode eval_mode; 7734b3e95d5SJeremy L Thompson 7742b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 7759e201c85SYohann if (eval_mode != CEED_EVAL_WEIGHT) { // Skip CEED_EVAL_WEIGHT 7764b3e95d5SJeremy L Thompson code << " const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n"; 7777d8d0e25Snbeams } 7787d8d0e25Snbeams } 7799e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 7804b3e95d5SJeremy L Thompson code << " CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n"; 7817d8d0e25Snbeams } 7827d8d0e25Snbeams 7839e201c85SYohann code << " const CeedInt dim = " << dim << ";\n"; 7849e201c85SYohann code << " const CeedInt Q_1d = " << Q_1d << ";\n"; 7857d8d0e25Snbeams 7864b3e95d5SJeremy L Thompson // Shared data 7874b3e95d5SJeremy L Thompson code << " extern __shared__ CeedScalar slice[];\n"; 7889e201c85SYohann code << " SharedData_Hip data;\n"; 7899e201c85SYohann code << " data.t_id_x = threadIdx.x;\n"; 7909e201c85SYohann code << " data.t_id_y = threadIdx.y;\n"; 7919e201c85SYohann code << " data.t_id_z = threadIdx.z;\n"; 7929e201c85SYohann code << " data.t_id = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n"; 7939e201c85SYohann code << " data.slice = slice + data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n"; 7947d8d0e25Snbeams 7957d8d0e25Snbeams // Initialize constants, and matrices B and G 7964b3e95d5SJeremy L Thompson code << "\n // Input field constants and basis data\n"; 7979e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 7985a5594ffSJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_input_fields[i], qf_input_fields[i], Q_1d, true, use_3d_slices)); 7997d8d0e25Snbeams } 8004b3e95d5SJeremy L Thompson code << "\n // Output field constants and basis data\n"; 8019e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 8025a5594ffSJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices)); 8034b3e95d5SJeremy L Thompson } 8047d8d0e25Snbeams 8054b3e95d5SJeremy L Thompson // Loop over all elements 8064b3e95d5SJeremy L Thompson code << "\n // Element loop\n"; 8077d8d0e25Snbeams code << " __syncthreads();\n"; 8089e201c85SYohann code << " for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n"; 8094b3e95d5SJeremy L Thompson 810e93651e5SJeremy L Thompson // -- Compute minimum buffer space needed 811e93651e5SJeremy L Thompson CeedInt max_rstr_buffer_size = 0; 812e93651e5SJeremy L Thompson 813e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 814e93651e5SJeremy L Thompson CeedInt num_comp, elem_size; 815e93651e5SJeremy L Thompson CeedElemRestriction elem_rstr; 816e93651e5SJeremy L Thompson 817e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 818e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 819e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 820e93651e5SJeremy L Thompson max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size); 821*681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 822e93651e5SJeremy L Thompson } 823e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 824e93651e5SJeremy L Thompson CeedInt num_comp, elem_size; 825e93651e5SJeremy L Thompson CeedElemRestriction elem_rstr; 826e93651e5SJeremy L Thompson 827e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 828e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 829e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 830e93651e5SJeremy L Thompson max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size); 831*681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 832e93651e5SJeremy L Thompson } 833e93651e5SJeremy L Thompson code << " // Scratch restriction buffer space\n"; 834e93651e5SJeremy L Thompson code << " CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n"; 835e93651e5SJeremy L Thompson 836e93651e5SJeremy L Thompson // -- Determine best input field processing order 837e93651e5SJeremy L Thompson CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX]; 838e93651e5SJeremy L Thompson 839e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 840e93651e5SJeremy L Thompson field_rstr_in_buffer[i] = -1; 841e93651e5SJeremy L Thompson input_field_order[i] = -1; 842e93651e5SJeremy L Thompson } 843e93651e5SJeremy L Thompson { 844e93651e5SJeremy L Thompson bool is_ordered[CEED_FIELD_MAX]; 845e93651e5SJeremy L Thompson CeedInt curr_index = 0; 846e93651e5SJeremy L Thompson 847e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false; 848e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 849e93651e5SJeremy L Thompson CeedVector vec_i; 850e93651e5SJeremy L Thompson CeedElemRestriction rstr_i; 851e93651e5SJeremy L Thompson 852e93651e5SJeremy L Thompson if (is_ordered[i]) continue; 853e93651e5SJeremy L Thompson field_rstr_in_buffer[i] = i; 854e93651e5SJeremy L Thompson is_ordered[i] = true; 855e93651e5SJeremy L Thompson input_field_order[curr_index] = i; 856e93651e5SJeremy L Thompson curr_index++; 857034f99fdSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i)); 858e93651e5SJeremy L Thompson if (vec_i == CEED_VECTOR_NONE) continue; // CEED_EVAL_WEIGHT 859e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i)); 860e93651e5SJeremy L Thompson for (CeedInt j = i + 1; j < num_input_fields; j++) { 861e93651e5SJeremy L Thompson CeedVector vec_j; 862e93651e5SJeremy L Thompson CeedElemRestriction rstr_j; 863e93651e5SJeremy L Thompson 864e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j)); 865e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j)); 866e93651e5SJeremy L Thompson if (rstr_i == rstr_j && vec_i == vec_j) { 867e93651e5SJeremy L Thompson field_rstr_in_buffer[j] = i; 868e93651e5SJeremy L Thompson is_ordered[j] = true; 869e93651e5SJeremy L Thompson input_field_order[curr_index] = j; 870e93651e5SJeremy L Thompson curr_index++; 871e93651e5SJeremy L Thompson } 872e93651e5SJeremy L Thompson } 873e93651e5SJeremy L Thompson } 874e93651e5SJeremy L Thompson } 875e93651e5SJeremy L Thompson 8764b3e95d5SJeremy L Thompson // -- Input restriction and basis 8774b3e95d5SJeremy L Thompson code << " // -- Input field restrictions and basis actions\n"; 8789e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 879e93651e5SJeremy L Thompson CeedInt f = input_field_order[i]; 880e93651e5SJeremy L Thompson 881e93651e5SJeremy L Thompson code << " // ---- Input field " << f << "\n"; 8827d8d0e25Snbeams 8834b3e95d5SJeremy L Thompson // ---- Restriction 884e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], Q_1d, 885e93651e5SJeremy L Thompson true, use_3d_slices)); 886b7453713SJeremy L Thompson 8874b3e95d5SJeremy L Thompson // ---- Basis action 888e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, use_3d_slices)); 8897d8d0e25Snbeams } 8907d8d0e25Snbeams 8914b3e95d5SJeremy L Thompson // -- Q function 8924b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelQFunction_Hip_gen(code, data, dim, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, 8934b3e95d5SJeremy L Thompson op_output_fields, qf_output_fields, qfunction_name, Q_1d, use_3d_slices)); 8947d8d0e25Snbeams 8954b3e95d5SJeremy L Thompson // -- Output basis and restriction 8964b3e95d5SJeremy L Thompson code << "\n // -- Output field basis action and restrictions\n"; 8979e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 8984b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 899b7453713SJeremy L Thompson 9004b3e95d5SJeremy L Thompson // ---- Basis action 9014b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices)); 9027d8d0e25Snbeams 9034b3e95d5SJeremy L Thompson // ---- Restriction 9044b3e95d5SJeremy L Thompson CeedCallBackend( 905e93651e5SJeremy L Thompson CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices)); 9067d8d0e25Snbeams } 9077d8d0e25Snbeams 9084b3e95d5SJeremy L Thompson // Close loop and function 9097d8d0e25Snbeams code << " }\n"; 9107d8d0e25Snbeams code << "}\n"; 9117d8d0e25Snbeams code << "// -----------------------------------------------------------------------------\n\n"; 9127d8d0e25Snbeams 9137d8d0e25Snbeams // View kernel for debugging 91423d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated Operator Kernels:\n"); 9153f21f6b1SJeremy L Thompson CeedDebug(ceed, code.str().c_str()); 9167d8d0e25Snbeams 917539ec17dSJeremy L Thompson CeedInt block_sizes[3] = {0, 0, 0}; 9189e201c85SYohann CeedInt num_elem; 919b7453713SJeremy L Thompson 9202b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 9212b730f8bSJeremy L Thompson CeedCallBackend(BlockGridCalculate_Hip_gen(dim, num_elem, data->max_P_1d, Q_1d, block_sizes)); 922eb7e6cafSJeremy L Thompson CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE", 9232b730f8bSJeremy L Thompson block_sizes[0] * block_sizes[1] * block_sizes[2])); 924eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op)); 9257d8d0e25Snbeams 9262b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorSetSetupDone(op)); 927e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 9287d8d0e25Snbeams } 9292a86cc9dSSebastian Grimberg 9307d8d0e25Snbeams //------------------------------------------------------------------------------ 931