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 struct FieldReuse_Hip { 25 CeedInt index; 26 bool is_input; 27 CeedEvalMode eval_mode; 28 }; 29 30 //------------------------------------------------------------------------------ 31 // Calculate the block size used for launching the operator kernel 32 //------------------------------------------------------------------------------ 33 extern "C" int BlockGridCalculate_Hip_gen(const CeedInt dim, const CeedInt num_elem, const CeedInt P_1d, const CeedInt Q_1d, CeedInt *block_sizes) { 34 const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 35 if (dim == 1) { 36 CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64; 37 38 elems_per_block = elems_per_block > 0 ? elems_per_block : 1; 39 block_sizes[0] = thread_1d; 40 block_sizes[1] = 1; 41 block_sizes[2] = elems_per_block; 42 } else if (dim == 2) { 43 const CeedInt elems_per_block = thread_1d < 4 ? 16 : 2; 44 45 block_sizes[0] = thread_1d; 46 block_sizes[1] = thread_1d; 47 block_sizes[2] = elems_per_block; 48 } else if (dim == 3) { 49 const CeedInt elems_per_block = thread_1d < 6 ? 4 : (thread_1d < 8 ? 2 : 1); 50 51 block_sizes[0] = thread_1d; 52 block_sizes[1] = thread_1d; 53 block_sizes[2] = elems_per_block; 54 } 55 return CEED_ERROR_SUCCESS; 56 } 57 58 //------------------------------------------------------------------------------ 59 // Determine type of operator 60 //------------------------------------------------------------------------------ 61 static int CeedOperatorBuildKernelData_Hip_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields, 62 CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields, 63 CeedQFunctionField *qf_output_fields, CeedInt *max_P, CeedInt *max_P_1d, CeedInt *Q, CeedInt *Q_1d, 64 CeedInt *max_dim, bool *is_all_tensor, bool *use_3d_slices) { 65 // Check if all are tensor 66 *is_all_tensor = true; 67 for (CeedInt i = 0; i < num_input_fields; i++) { 68 CeedBasis basis; 69 70 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 71 if (basis != CEED_BASIS_NONE) { 72 bool is_field_tensor; 73 74 CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 75 *is_all_tensor = *is_all_tensor && is_field_tensor; 76 } 77 CeedCallBackend(CeedBasisDestroy(&basis)); 78 } 79 for (CeedInt i = 0; i < num_output_fields; i++) { 80 CeedBasis basis; 81 82 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 83 if (basis != CEED_BASIS_NONE) { 84 bool is_field_tensor; 85 86 CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 87 *is_all_tensor = *is_all_tensor && is_field_tensor; 88 } 89 CeedCallBackend(CeedBasisDestroy(&basis)); 90 } 91 92 // Find max_P, max_P_1d, Q, and Q_1d 93 bool is_all_3d = true; 94 95 *max_P = 0; 96 *max_P_1d = 0; 97 *Q = 0; 98 *Q_1d = 0; 99 for (CeedInt i = 0; i < num_input_fields; i++) { 100 CeedBasis basis; 101 102 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 103 if (basis != CEED_BASIS_NONE) { 104 bool is_field_tensor; 105 CeedInt field_dim = 0, field_P = 0, field_P_1d = 0, field_Q = 0, field_Q_1d = 0; 106 107 // Check if 3D 108 CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 109 is_all_3d = is_all_3d && (field_dim == 3); 110 *max_dim = CeedIntMax(*max_dim, field_dim); 111 112 // Collect P, P_1d, Q, and Q_1d 113 CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P)); 114 *max_P = CeedIntMax(*max_P, field_P); 115 CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 116 if (is_field_tensor) { 117 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 118 *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 119 } 120 CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q)); 121 CeedCheck(*Q == 0 || field_Q == *Q, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 122 *Q = field_Q; 123 if (is_field_tensor) { 124 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 125 CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 126 *Q_1d = field_Q_1d; 127 } 128 } 129 CeedCallBackend(CeedBasisDestroy(&basis)); 130 } 131 for (CeedInt i = 0; i < num_output_fields; i++) { 132 CeedBasis basis; 133 134 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 135 if (basis != CEED_BASIS_NONE) { 136 bool is_field_tensor; 137 CeedInt field_dim = 0, field_P = 0, field_P_1d = 0, field_Q = 0, field_Q_1d = 0; 138 139 // Check if 3D 140 CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 141 is_all_3d = is_all_3d && (field_dim == 3); 142 *max_dim = CeedIntMax(*max_dim, field_dim); 143 144 // Collect P, P_1d, Q, and Q_1d 145 CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P)); 146 *max_P = CeedIntMax(*max_P, field_P); 147 CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 148 if (is_field_tensor) { 149 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 150 *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 151 } 152 CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q)); 153 CeedCheck(*Q == 0 || field_Q == *Q, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 154 *Q = field_Q; 155 if (is_field_tensor) { 156 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 157 CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 158 *Q_1d = field_Q_1d; 159 } 160 } 161 CeedCallBackend(CeedBasisDestroy(&basis)); 162 } 163 164 // Only use 3D collocated gradient parallelization strategy when gradient is computed 165 *use_3d_slices = false; 166 if (is_all_3d && *is_all_tensor) { 167 bool was_grad_found = false; 168 169 for (CeedInt i = 0; i < num_input_fields; i++) { 170 CeedEvalMode eval_mode; 171 172 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 173 if (eval_mode == CEED_EVAL_GRAD) { 174 CeedBasis_Hip_shared *basis_data; 175 CeedBasis basis; 176 177 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 178 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 179 *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 180 was_grad_found = true; 181 CeedCallBackend(CeedBasisDestroy(&basis)); 182 } 183 } 184 for (CeedInt i = 0; i < num_output_fields; i++) { 185 CeedEvalMode eval_mode; 186 187 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 188 if (eval_mode == CEED_EVAL_GRAD) { 189 CeedBasis_Hip_shared *basis_data; 190 CeedBasis basis; 191 192 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 193 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 194 *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 195 was_grad_found = true; 196 CeedCallBackend(CeedBasisDestroy(&basis)); 197 } 198 } 199 } 200 return CEED_ERROR_SUCCESS; 201 } 202 203 //------------------------------------------------------------------------------ 204 // Setup fields 205 //------------------------------------------------------------------------------ 206 static int CeedOperatorBuildKernelFieldData_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedOperatorField op_field, 207 CeedQFunctionField qf_field, FieldReuse_Hip field_reuse, CeedInt max_dim, CeedInt Q, CeedInt Q_1d, 208 bool is_input, bool is_all_tensor, bool is_at_points, bool use_3d_slices) { 209 bool is_tensor = true; 210 CeedBasis basis; 211 CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 212 if (basis != CEED_BASIS_NONE) CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor)); 213 214 const char *field_name; 215 std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 216 std::string P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q"; 217 std::string option_name = (is_input ? "inputs" : "outputs"); 218 CeedEvalMode eval_mode = CEED_EVAL_NONE; 219 CeedInt elem_size = 0, num_comp = 0, dim = max_dim, P_1d = 0; 220 CeedElemRestriction elem_rstr; 221 CeedBasis_Hip_shared *basis_data; 222 223 // Field reuse info 224 bool use_previous_field = field_reuse.index != -1; 225 226 CeedCallBackend(CeedOperatorFieldGetName(op_field, &field_name)); 227 code << " // -- " << (is_input ? "Input" : "Output") << " field " << i << ": " << field_name << "\n"; 228 229 // Get field data 230 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 231 if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 232 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 233 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 234 } 235 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 236 if (basis != CEED_BASIS_NONE) { 237 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 238 CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 239 if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 240 else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d)); 241 } 242 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 243 244 // Set field constants 245 code << " const CeedInt dim" << var_suffix << " = " << dim << ";\n"; 246 if (is_tensor && !is_all_tensor) { 247 CeedInt P = 0; 248 249 CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 250 code << " const CeedInt P" << var_suffix << " = " << (basis == CEED_BASIS_NONE ? Q : P) << ";\n"; 251 } 252 code << " const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n"; 253 if (eval_mode != CEED_EVAL_WEIGHT) { 254 code << " const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n"; 255 } 256 257 // Load basis data 258 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 259 switch (eval_mode) { 260 case CEED_EVAL_NONE: 261 break; 262 case CEED_EVAL_INTERP: 263 if (is_at_points) { 264 // AtPoints 265 if (!basis_data->d_chebyshev_interp_1d) { 266 CeedSize interp_bytes; 267 CeedScalar *chebyshev_interp_1d; 268 269 interp_bytes = P_1d * Q_1d * sizeof(CeedScalar); 270 CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 271 CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 272 CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes)); 273 CeedCallHip(CeedBasisReturnCeed(basis), 274 hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice)); 275 CeedCallBackend(CeedFree(&chebyshev_interp_1d)); 276 } 277 if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d; 278 else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d; 279 } else { 280 // Standard quadrature 281 if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 282 else data->B.outputs[i] = basis_data->d_interp_1d; 283 } 284 if (use_previous_field) { 285 std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index)); 286 287 code << " CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n"; 288 } else { 289 code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n"; 290 code << " LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 291 } 292 break; 293 case CEED_EVAL_GRAD: 294 if (is_at_points) { 295 // AtPoints 296 if (!basis_data->d_chebyshev_interp_1d) { 297 CeedSize interp_bytes; 298 CeedScalar *chebyshev_interp_1d; 299 300 interp_bytes = P_1d * Q_1d * sizeof(CeedScalar); 301 CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 302 CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 303 CeedCallHip(CeedBasisReturnCeed(basis), hipMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes)); 304 CeedCallHip(CeedBasisReturnCeed(basis), 305 hipMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, hipMemcpyHostToDevice)); 306 CeedCallBackend(CeedFree(&chebyshev_interp_1d)); 307 } 308 if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d; 309 else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d; 310 } else { 311 // Standard quadrature 312 if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 313 else data->B.outputs[i] = basis_data->d_interp_1d; 314 } 315 if (is_tensor) { 316 if (use_previous_field) { 317 std::string reuse_var = "s_B" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index)); 318 319 code << " CeedScalar *s_B" << var_suffix << " = " << reuse_var << ";\n"; 320 } else { 321 code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n"; 322 code << " LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 323 } 324 } 325 if (is_at_points) break; // No G mat for AtPoints 326 if (use_3d_slices) { 327 if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d; 328 else data->G.outputs[i] = basis_data->d_collo_grad_1d; 329 if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) { 330 std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index)); 331 332 code << " CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n"; 333 } else { 334 code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n"; 335 code << " LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 336 } 337 } else { 338 bool has_collo_grad = basis_data->d_collo_grad_1d; 339 340 if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 341 else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 342 if (has_collo_grad) { 343 if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) { 344 std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index)); 345 346 code << " CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n"; 347 } else { 348 code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n"; 349 code << " LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 350 } 351 } else { 352 if (use_previous_field && field_reuse.eval_mode == CEED_EVAL_GRAD) { 353 std::string reuse_var = "s_G" + ((field_reuse.is_input ? "_in_" : "_out_") + std::to_string(field_reuse.index)); 354 355 code << " CeedScalar *s_G" << var_suffix << " = " << reuse_var << ";\n"; 356 } else { 357 code << " __shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim") 358 << (is_tensor ? "" : var_suffix) << "];\n"; 359 code << " LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << (is_tensor ? "" : var_suffix) << ">(data, G." 360 << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 361 } 362 } 363 } 364 break; 365 case CEED_EVAL_WEIGHT: 366 break; // No action 367 // LCOV_EXCL_START 368 case CEED_EVAL_DIV: 369 case CEED_EVAL_CURL: 370 break; // TODO: Not implemented 371 // LCOV_EXCL_STOP 372 } 373 CeedCallBackend(CeedBasisDestroy(&basis)); 374 return CEED_ERROR_SUCCESS; 375 } 376 377 //------------------------------------------------------------------------------ 378 // Restriction 379 //------------------------------------------------------------------------------ 380 static int CeedOperatorBuildKernelRestriction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedInt field_input_buffer[], 381 CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt max_dim, CeedInt Q_1d, 382 bool is_input, bool is_all_tensor, bool is_at_points, bool use_3d_slices) { 383 std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 384 std::string P_name = (is_all_tensor ? "P_1d" : "P") + var_suffix; 385 CeedEvalMode eval_mode = CEED_EVAL_NONE; 386 CeedInt elem_size = 0, num_comp = 0; 387 CeedSize l_size; 388 CeedRestrictionType rstr_type = CEED_RESTRICTION_STANDARD; 389 CeedElemRestriction_Hip *rstr_data; 390 CeedElemRestriction elem_rstr; 391 392 // Get field data 393 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 394 if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 395 CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type)); 396 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 397 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 398 CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 399 } 400 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 401 402 // Restriction 403 if (is_input) { 404 // Input 405 if (field_input_buffer[i] != i) { 406 std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]); 407 408 // Restriction was already done for previous input 409 code << " CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n"; 410 } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) { 411 if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) { 412 // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated 413 code << " CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n"; 414 } else if (rstr_type != CEED_RESTRICTION_POINTS) { 415 // Otherwise we're using the scratch space 416 code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 417 } 418 switch (rstr_type) { 419 case CEED_RESTRICTION_STANDARD: { 420 CeedInt comp_stride; 421 422 CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 423 code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 424 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 425 code << " // CompStride: " << comp_stride << "\n"; 426 data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 427 code << " ReadLVecStandard" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name 428 << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n"; 429 break; 430 } 431 case CEED_RESTRICTION_STRIDED: { 432 bool has_backend_strides; 433 CeedInt num_elem; 434 435 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 436 CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 437 CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 438 439 if (!has_backend_strides) { 440 CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 441 } 442 code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 443 code << " ReadLVecStrided" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", " 444 << strides[1] << ", " << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n"; 445 break; 446 } 447 case CEED_RESTRICTION_POINTS: { 448 CeedInt comp_stride; 449 450 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 451 code << " const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n"; 452 data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 453 break; 454 } 455 // LCOV_EXCL_START 456 case CEED_RESTRICTION_ORIENTED: 457 case CEED_RESTRICTION_CURL_ORIENTED: 458 break; // TODO: Not implemented 459 // LCOV_EXCL_STOP 460 } 461 } 462 } else { 463 // Output 464 switch (rstr_type) { 465 case CEED_RESTRICTION_STANDARD: { 466 CeedInt comp_stride; 467 468 CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 469 code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 470 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 471 code << " // CompStride: " << comp_stride << "\n"; 472 data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets; 473 code << " WriteLVecStandard" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name 474 << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n"; 475 break; 476 } 477 case CEED_RESTRICTION_STRIDED: { 478 bool has_backend_strides; 479 CeedInt num_elem; 480 481 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 482 CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 483 CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 484 485 if (!has_backend_strides) { 486 CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 487 } 488 code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 489 code << " WriteLVecStrided" << (is_all_tensor ? max_dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", " 490 << strides[1] << ", " << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n"; 491 break; 492 } 493 case CEED_RESTRICTION_POINTS: 494 data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets; 495 break; 496 // LCOV_EXCL_START 497 case CEED_RESTRICTION_ORIENTED: 498 case CEED_RESTRICTION_CURL_ORIENTED: 499 break; // TODO: Not implemented 500 // LCOV_EXCL_STOP 501 } 502 } 503 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 504 return CEED_ERROR_SUCCESS; 505 } 506 507 //------------------------------------------------------------------------------ 508 // Basis 509 //------------------------------------------------------------------------------ 510 static int CeedOperatorBuildKernelBasis_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt i, CeedOperatorField op_field, 511 CeedQFunctionField qf_field, CeedInt max_dim, CeedInt Q_1d, bool is_input, bool is_all_tensor, 512 bool is_at_points, bool use_3d_slices) { 513 bool is_tensor = true; 514 CeedBasis basis; 515 CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 516 CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor)); 517 518 std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 519 std::string P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q"; 520 CeedEvalMode eval_mode = CEED_EVAL_NONE; 521 CeedInt dim = max_dim, elem_size = 0, num_comp = 0, P_1d = 0; 522 CeedElemRestriction elem_rstr; 523 524 // Get field data 525 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 526 if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 527 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 528 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 529 } 530 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 531 if (basis != CEED_BASIS_NONE) { 532 CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 533 if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 534 else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d)); 535 } 536 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 537 538 // Basis 539 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 540 if (is_input) { 541 switch (eval_mode) { 542 case CEED_EVAL_NONE: 543 if (!use_3d_slices && !is_at_points) { 544 code << " CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n"; 545 } 546 break; 547 case CEED_EVAL_INTERP: 548 if (is_at_points) { 549 std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d"; 550 551 code << " CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n"; 552 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix 553 << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n"; 554 } else { 555 std::string function_name = is_tensor 556 ? ((dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened")) 557 : "InterpNonTensor"; 558 std::string op_t_1d_name = (is_all_tensor || !is_tensor) ? "OP_T_1D" : (P_1d > Q_1d ? P_name : Q_name); 559 560 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_all_tensor && (dim >= 3) ? Q_name : "1") << "];\n"; 561 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << op_t_1d_name << ">(data, r_e" 562 << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n"; 563 } 564 break; 565 case CEED_EVAL_GRAD: 566 if (is_at_points) { 567 std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d"; 568 569 code << " CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n"; 570 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix 571 << ", s_B" << var_suffix << ", r_c" << var_suffix << ");\n"; 572 } else if (use_3d_slices) { 573 std::string function_name = (dim > 1 ? "InterpTensor" : "Interp") + std::to_string(dim) + "d"; 574 575 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 576 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_e" << var_suffix 577 << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n"; 578 } else if (is_tensor) { 579 bool is_collocated = dim == 3 && Q_1d >= P_1d; 580 std::string function_name = (dim == 1 ? "Grad" : (is_collocated ? "GradTensorCollocated" : "GradTensor")) + std::to_string(dim) + "d" + 581 (is_all_tensor ? "" : "Flattened"); 582 std::string op_t_1d_name = is_all_tensor ? "OP_T_1D" : (P_1d > Q_1d ? P_name : Q_name); 583 584 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "*" 585 << (is_all_tensor && dim >= 3 ? Q_name : "1") << "];\n"; 586 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << op_t_1d_name << ">(data, r_e" 587 << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n"; 588 } else { 589 std::string function_name = "GradNonTensor"; 590 591 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n"; 592 code << " " << function_name << "<num_comp" << var_suffix << ", dim" << var_suffix << ", " << P_name << ", " << Q_name 593 << ", OP_T_1D>(data, r_e" << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n"; 594 } 595 break; 596 case CEED_EVAL_WEIGHT: { 597 if (is_at_points) { 598 code << " // Nothing to do AtPoints\n"; 599 } else { 600 CeedBasis_Hip_shared *basis_data; 601 std::string function_name = is_tensor 602 ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened")) 603 : "WeightNonTensor"; 604 605 code << " CeedScalar r_q" << var_suffix << "[" << (is_all_tensor && (dim >= 3) ? Q_name : "1") << "];\n"; 606 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 607 data->W = basis_data->d_q_weight_1d; 608 code << " " << function_name << "<" << P_name << ", " << Q_name << ">(data, W, r_q" << var_suffix << ");\n"; 609 } 610 break; 611 } 612 // LCOV_EXCL_START 613 case CEED_EVAL_DIV: 614 case CEED_EVAL_CURL: 615 break; // TODO: Not implemented 616 // LCOV_EXCL_STOP 617 } 618 } else { 619 switch (eval_mode) { 620 case CEED_EVAL_NONE: 621 code << " CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n"; 622 break; // No action 623 case CEED_EVAL_INTERP: 624 code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 625 if (is_at_points) { 626 std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d"; 627 628 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_c" << var_suffix 629 << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 630 } else { 631 std::string function_name = 632 is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened")) 633 : "InterpTransposeNonTensor"; 634 std::string op_t_1d_name = (is_all_tensor || !is_tensor) ? "OP_T_1D" : (P_1d > Q_1d ? P_name : Q_name); 635 636 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << op_t_1d_name << ">(data, r_q" 637 << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 638 } 639 break; 640 case CEED_EVAL_GRAD: 641 code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 642 if (is_at_points) { 643 std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d"; 644 645 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_c" << var_suffix 646 << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 647 } else if (use_3d_slices) { 648 std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d"; 649 650 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", OP_T_1D>(data, r_q" << var_suffix 651 << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 652 } else if (is_tensor) { 653 bool is_collocated = dim == 3 && Q_1d >= P_1d; 654 std::string function_name = (dim == 1 ? "GradTranspose" : (is_collocated ? "GradTransposeTensorCollocated" : "GradTransposeTensor")) + 655 std::to_string(dim) + "d" + (is_all_tensor ? "" : "Flattened"); 656 std::string op_t_1d_name = is_all_tensor ? "OP_T_1D" : (P_1d > Q_1d ? P_name : Q_name); 657 658 code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ", " << op_t_1d_name << ">(data, r_q" 659 << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n"; 660 } else { 661 std::string function_name = "GradTransposeNonTensor"; 662 663 code << " " << function_name << "<num_comp" << var_suffix << ", dim" << var_suffix << ", " << P_name << ", " << Q_name 664 << ", OP_T_1D>(data, r_q" << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n"; 665 } 666 break; 667 // LCOV_EXCL_START 668 case CEED_EVAL_WEIGHT: 669 break; // Should not occur 670 case CEED_EVAL_DIV: 671 case CEED_EVAL_CURL: 672 break; // TODO: Not implemented 673 // LCOV_EXCL_STOP 674 } 675 } 676 CeedCallBackend(CeedBasisDestroy(&basis)); 677 return CEED_ERROR_SUCCESS; 678 } 679 680 //------------------------------------------------------------------------------ 681 // QFunction 682 //------------------------------------------------------------------------------ 683 static int CeedOperatorBuildKernelQFunction_Hip_gen(std::ostringstream &code, CeedOperator_Hip_gen *data, CeedInt max_dim, CeedInt max_num_points, 684 CeedInt num_input_fields, CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields, 685 CeedInt num_output_fields, CeedOperatorField *op_output_fields, 686 CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d, 687 bool is_all_tensor, bool is_at_points, bool use_3d_slices) { 688 std::string Q_name = is_all_tensor ? "Q_1d" : "Q"; 689 CeedEvalMode eval_mode = CEED_EVAL_NONE; 690 CeedElemRestriction elem_rstr; 691 692 // Setup output arrays 693 code << "\n // -- Output field setup\n"; 694 for (CeedInt i = 0; i < num_output_fields; i++) { 695 const char *field_name; 696 std::string var_suffix = "_out_" + std::to_string(i); 697 698 CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name)); 699 code << " // ---- Output field " << i << ": " << field_name << "\n"; 700 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 701 switch (eval_mode) { 702 case CEED_EVAL_NONE: 703 if (is_at_points) { 704 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n"; 705 } else { 706 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_all_tensor && (max_dim >= 3) ? Q_name : "1") 707 << "];\n"; 708 } 709 break; 710 case CEED_EVAL_INTERP: 711 if (is_at_points) { 712 // Accumulator for point data 713 code << " CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "];\n"; 714 code << " for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "; i++) {\n"; 715 code << " r_c" << var_suffix << "[i] = 0.0;\n"; 716 code << " }\n"; 717 } else { 718 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_all_tensor && (max_dim >= 3) ? Q_name : "1") 719 << "];\n"; 720 } 721 break; 722 case CEED_EVAL_GRAD: 723 if (is_at_points) { 724 // Accumulator for point data 725 code << " CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "*dim" << var_suffix 726 << "];\n"; 727 code << " for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (max_dim >= 3 ? Q_name : "1") << "; i++) {\n"; 728 code << " r_c" << var_suffix << "[i] = 0.0;\n"; 729 code << " }\n"; 730 } else if (use_3d_slices) { 731 // Accumulator for gradient slices 732 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 733 code << " for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n"; 734 code << " r_q" << var_suffix << "[i] = 0.0;\n"; 735 code << " }\n"; 736 } else { 737 code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "*" 738 << (is_all_tensor && (max_dim >= 3) ? Q_name : "1") << "];\n"; 739 } 740 break; 741 case CEED_EVAL_WEIGHT: 742 break; 743 // LCOV_EXCL_START 744 case CEED_EVAL_DIV: 745 case CEED_EVAL_CURL: 746 break; // TODO: Not implemented 747 // LCOV_EXCL_STOP 748 } 749 } 750 751 if (is_at_points) { 752 // We need to handle batches of points 753 code << "\n // Note: Using batches of points\n"; 754 code << " const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n"; 755 code << " #pragma unroll\n"; 756 code << " for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n"; 757 code << " const CeedInt p = i % max_num_points;\n\n"; 758 759 code << " // -- Coordinates\n"; 760 code << " CeedScalar r_x[max_dim];\n"; 761 code << " ReadPoint<max_dim, coords_comp_stride, max_num_points>(data, elem, p, max_num_points, points.indices, points.coords, r_x);\n\n"; 762 763 code << " // -- Input fields\n"; 764 for (CeedInt i = 0; i < num_input_fields; i++) { 765 const char *field_name; 766 std::string var_suffix = "_in_" + std::to_string(i); 767 std::string P_name = "P_1d" + var_suffix; 768 769 CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name)); 770 code << " // ---- Input field " << i << ": " << field_name << "\n"; 771 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 772 // Basis action 773 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 774 switch (eval_mode) { 775 case CEED_EVAL_NONE: 776 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 777 code << " ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix 778 << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n"; 779 break; 780 case CEED_EVAL_INTERP: 781 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 782 code << " InterpAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name 783 << ">(data, i, r_c" << var_suffix << ", r_x, r_s" << var_suffix << ");\n"; 784 break; 785 case CEED_EVAL_GRAD: 786 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n"; 787 code << " GradAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name 788 << ">(data, i, r_c" << var_suffix << ", r_x, r_s" << var_suffix << ");\n"; 789 break; 790 case CEED_EVAL_WEIGHT: 791 code << " CeedScalar r_s" << var_suffix << "[1];\n"; 792 code << " r_s" << var_suffix << "[0] = 1.0;\n"; 793 break; 794 // LCOV_EXCL_START 795 case CEED_EVAL_DIV: 796 case CEED_EVAL_CURL: 797 break; // TODO: Not implemented 798 // LCOV_EXCL_STOP 799 } 800 } 801 code << "\n // -- Output fields\n"; 802 for (CeedInt i = 0; i < num_output_fields; i++) { 803 const char *field_name; 804 std::string var_suffix = "_out_" + std::to_string(i); 805 806 CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name)); 807 code << " // ---- Output field " << i << ": " << field_name << "\n"; 808 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 809 // Basis action 810 switch (eval_mode) { 811 case CEED_EVAL_NONE: 812 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 813 break; 814 case CEED_EVAL_INTERP: 815 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 816 break; 817 case CEED_EVAL_GRAD: 818 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n"; 819 break; 820 // LCOV_EXCL_START 821 case CEED_EVAL_WEIGHT: 822 break; // Should not occur 823 case CEED_EVAL_DIV: 824 case CEED_EVAL_CURL: 825 break; // TODO: Not implemented 826 // LCOV_EXCL_STOP 827 } 828 } 829 830 } else if (use_3d_slices) { 831 // We treat quadrature points per slice in 3d to save registers 832 code << "\n // Note: Using planes of 3D elements\n"; 833 code << " #pragma unroll\n"; 834 code << " for (CeedInt q = 0; q < " << Q_name << "; q++) {\n"; 835 code << " // -- Input fields\n"; 836 for (CeedInt i = 0; i < num_input_fields; i++) { 837 const char *field_name; 838 std::string var_suffix = "_in_" + std::to_string(i); 839 840 CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name)); 841 code << " // ---- Input field " << i << ": " << field_name << "\n"; 842 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 843 // Basis action 844 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 845 switch (eval_mode) { 846 case CEED_EVAL_NONE: 847 bool is_strided; 848 849 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 850 851 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 852 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 853 if (is_strided) { 854 bool has_backend_strides; 855 CeedInt num_elem, elem_size; 856 857 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 858 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 859 CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 860 CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 861 862 if (!has_backend_strides) { 863 CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 864 } 865 code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 866 code << " ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", " << strides[0] << ", " << strides[1] << ", " 867 << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n"; 868 } else { 869 CeedSize l_size = 0; 870 CeedInt comp_stride; 871 CeedElemRestriction_Hip *rstr_data; 872 873 CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 874 code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 875 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 876 code << " // CompStride: " << comp_stride << "\n"; 877 CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 878 data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 879 code << " ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix 880 << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n"; 881 } 882 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 883 break; 884 case CEED_EVAL_INTERP: 885 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 886 code << " for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n"; 887 code << " r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n"; 888 code << " }\n"; 889 break; 890 case CEED_EVAL_GRAD: 891 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n"; 892 code << " GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ", OP_T_1D>(data, q, r_q" << var_suffix << ", s_G" 893 << var_suffix << ", r_s" << var_suffix << ");\n"; 894 break; 895 case CEED_EVAL_WEIGHT: 896 code << " CeedScalar r_s" << var_suffix << "[1];\n"; 897 code << " r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n"; 898 break; 899 // LCOV_EXCL_START 900 case CEED_EVAL_DIV: 901 case CEED_EVAL_CURL: 902 break; // TODO: Not implemented 903 // LCOV_EXCL_STOP 904 } 905 } 906 code << "\n // -- Output fields\n"; 907 for (CeedInt i = 0; i < num_output_fields; i++) { 908 const char *field_name; 909 std::string var_suffix = "_out_" + std::to_string(i); 910 911 CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name)); 912 code << " // ---- Output field " << i << ": " << field_name << "\n"; 913 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 914 // Basis action 915 switch (eval_mode) { 916 case CEED_EVAL_NONE: 917 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 918 break; 919 case CEED_EVAL_INTERP: 920 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 921 break; 922 case CEED_EVAL_GRAD: 923 code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim" << var_suffix << "];\n"; 924 break; 925 // LCOV_EXCL_START 926 case CEED_EVAL_WEIGHT: 927 break; // Should not occur 928 case CEED_EVAL_DIV: 929 case CEED_EVAL_CURL: 930 break; // TODO: Not implemented 931 // LCOV_EXCL_STOP 932 } 933 } 934 } else { 935 code << "\n // Note: Using full elements\n"; 936 code << " {\n"; 937 code << " // -- Input fields\n"; 938 for (CeedInt i = 0; i < num_input_fields; i++) { 939 const char *field_name; 940 941 CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name)); 942 code << " // ---- Input field " << i << ": " << field_name << "\n"; 943 code << " CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n"; 944 } 945 code << " // -- Output fields\n"; 946 for (CeedInt i = 0; i < num_output_fields; i++) { 947 const char *field_name; 948 949 CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name)); 950 code << " // ---- Output field " << i << ": " << field_name << "\n"; 951 code << " CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n"; 952 } 953 } 954 955 // Input and output buffers 956 code << "\n // -- QFunction inputs and outputs\n"; 957 code << " // ---- Inputs\n"; 958 code << " CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n"; 959 for (CeedInt i = 0; i < num_input_fields; i++) { 960 const char *field_name; 961 962 CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[i], &field_name)); 963 code << " // ------ Input field " << i << ": " << field_name << "\n"; 964 code << " inputs[" << i << "] = r_s_in_" << i << ";\n"; 965 } 966 code << " // ---- Outputs\n"; 967 code << " CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n"; 968 for (CeedInt i = 0; i < num_output_fields; i++) { 969 const char *field_name; 970 971 CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name)); 972 code << " // ------ Output field " << i << ": " << field_name << "\n"; 973 code << " outputs[" << i << "] = r_s_out_" << i << ";\n"; 974 } 975 976 // Apply QFunction 977 code << "\n // -- Apply QFunction\n"; 978 code << " " << qfunction_name << "(ctx, "; 979 if (max_dim != 3 || is_at_points || use_3d_slices || !is_all_tensor) { 980 code << "1"; 981 } else { 982 code << Q_name; 983 } 984 code << ", inputs, outputs);\n"; 985 986 if (is_at_points) { 987 // Map back to coefficients 988 code << "\n // -- Output fields\n"; 989 for (CeedInt i = 0; i < num_output_fields; i++) { 990 const char *field_name; 991 std::string var_suffix = "_out_" + std::to_string(i); 992 std::string P_name = "P_1d" + var_suffix; 993 994 CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name)); 995 code << " // ---- Output field " << i << ": " << field_name << "\n"; 996 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 997 // Basis action 998 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 999 switch (eval_mode) { 1000 case CEED_EVAL_NONE: { 1001 CeedInt comp_stride; 1002 CeedElemRestriction elem_rstr; 1003 1004 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 1005 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 1006 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 1007 code << " const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n"; 1008 code << " WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix 1009 << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]" 1010 << ", r_s" << var_suffix << ", d" << var_suffix << ");\n"; 1011 break; 1012 } 1013 case CEED_EVAL_INTERP: 1014 code << " if (i >= points.num_per_elem[elem]) {\n"; 1015 code << " for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n"; 1016 code << " }\n"; 1017 code << " InterpTransposeAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name 1018 << ">(data, i, r_s" << var_suffix << ", r_x, r_c" << var_suffix << ");\n"; 1019 break; 1020 case CEED_EVAL_GRAD: 1021 code << " if (i >= points.num_per_elem[elem]) {\n"; 1022 code << " for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n"; 1023 code << " }\n"; 1024 code << " GradTransposeAtPoints" << max_dim << "d<num_comp" << var_suffix << ", max_num_points, " << P_name << ", " << Q_name 1025 << ">(data, i, r_s" << var_suffix << ", r_x, r_c" << var_suffix << ");\n"; 1026 break; 1027 // LCOV_EXCL_START 1028 case CEED_EVAL_WEIGHT: 1029 break; // Should not occur 1030 case CEED_EVAL_DIV: 1031 case CEED_EVAL_CURL: 1032 break; // TODO: Not implemented 1033 // LCOV_EXCL_STOP 1034 } 1035 } 1036 } else if (use_3d_slices) { 1037 // Copy or apply transpose grad, if needed 1038 code << "\n // -- Output fields\n"; 1039 for (CeedInt i = 0; i < num_output_fields; i++) { 1040 const char *field_name; 1041 std::string var_suffix = "_out_" + std::to_string(i); 1042 std::string P_name = "P_1d" + var_suffix; 1043 1044 CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name)); 1045 code << " // ---- Output field " << i << ": " << field_name << "\n"; 1046 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1047 // Basis action 1048 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 1049 switch (eval_mode) { 1050 case CEED_EVAL_NONE: 1051 code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 1052 code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 1053 code << " }\n"; 1054 break; 1055 case CEED_EVAL_INTERP: 1056 code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 1057 code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 1058 code << " }\n"; 1059 break; 1060 case CEED_EVAL_GRAD: 1061 code << " GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ", OP_T_1D>(data, q, r_s" << var_suffix << ", s_G" 1062 << var_suffix << ", r_q" << var_suffix << ");\n"; 1063 break; 1064 // LCOV_EXCL_START 1065 case CEED_EVAL_WEIGHT: 1066 break; // Should not occur 1067 case CEED_EVAL_DIV: 1068 case CEED_EVAL_CURL: 1069 break; // TODO: Not implemented 1070 // LCOV_EXCL_STOP 1071 } 1072 } 1073 } 1074 code << " }\n"; 1075 return CEED_ERROR_SUCCESS; 1076 } 1077 1078 //------------------------------------------------------------------------------ 1079 // Build single operator kernel 1080 //------------------------------------------------------------------------------ 1081 extern "C" int CeedOperatorBuildKernel_Hip_gen(CeedOperator op, bool *is_good_build) { 1082 bool is_all_tensor = true, is_all_nontensor = true, is_at_points = false, use_3d_slices = false; 1083 Ceed ceed; 1084 CeedInt Q, Q_1d, num_input_fields, num_output_fields, max_dim = 1, max_num_points = 0, coords_comp_stride = 0; 1085 CeedQFunctionField *qf_input_fields, *qf_output_fields; 1086 CeedQFunction_Hip_gen *qf_data; 1087 CeedQFunction qf; 1088 CeedOperatorField *op_input_fields, *op_output_fields; 1089 CeedOperator_Hip_gen *data; 1090 std::ostringstream code; 1091 1092 CeedCallBackend(CeedOperatorGetData(op, &data)); 1093 { 1094 bool is_setup_done; 1095 1096 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 1097 if (is_setup_done) { 1098 *is_good_build = !data->use_fallback; 1099 return CEED_ERROR_SUCCESS; 1100 } 1101 } 1102 1103 // Check field compatibility 1104 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 1105 { 1106 bool has_shared_bases = true; 1107 1108 for (CeedInt i = 0; i < num_input_fields; i++) { 1109 CeedBasis basis; 1110 1111 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 1112 if (basis != CEED_BASIS_NONE) { 1113 bool is_tensor = true; 1114 const char *resource; 1115 char *resource_root; 1116 Ceed basis_ceed; 1117 1118 CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor)); 1119 is_all_tensor = is_all_tensor && is_tensor; 1120 is_all_nontensor = is_all_nontensor && !is_tensor; 1121 CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed)); 1122 CeedCallBackend(CeedGetResource(basis_ceed, &resource)); 1123 CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root)); 1124 has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/hip/shared"); 1125 CeedCallBackend(CeedFree(&resource_root)); 1126 CeedCallBackend(CeedDestroy(&basis_ceed)); 1127 } 1128 CeedCallBackend(CeedBasisDestroy(&basis)); 1129 } 1130 1131 for (CeedInt i = 0; i < num_output_fields; i++) { 1132 CeedBasis basis; 1133 1134 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 1135 if (basis != CEED_BASIS_NONE) { 1136 bool is_tensor = true; 1137 const char *resource; 1138 char *resource_root; 1139 Ceed basis_ceed; 1140 1141 CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor)); 1142 is_all_tensor = is_all_tensor && is_tensor; 1143 is_all_nontensor = is_all_nontensor && !is_tensor; 1144 1145 CeedCallBackend(CeedBasisGetCeed(basis, &basis_ceed)); 1146 CeedCallBackend(CeedGetResource(basis_ceed, &resource)); 1147 CeedCallBackend(CeedGetResourceRoot(basis_ceed, resource, ":", &resource_root)); 1148 has_shared_bases = has_shared_bases && !strcmp(resource_root, "/gpu/hip/shared"); 1149 CeedCallBackend(CeedFree(&resource_root)); 1150 CeedCallBackend(CeedDestroy(&basis_ceed)); 1151 } 1152 CeedCallBackend(CeedBasisDestroy(&basis)); 1153 } 1154 // -- Fallback to ref if not all bases are shared 1155 if (!has_shared_bases) { 1156 *is_good_build = false; 1157 return CEED_ERROR_SUCCESS; 1158 } 1159 } 1160 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1161 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 1162 CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 1163 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 1164 1165 // Get operator data 1166 CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points)); 1167 { 1168 CeedInt max_P, max_P_1d; 1169 1170 CeedCallBackend(CeedOperatorBuildKernelData_Hip_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields, 1171 qf_output_fields, &max_P, &max_P_1d, &Q, &Q_1d, &max_dim, &is_all_tensor, &use_3d_slices)); 1172 data->max_P_1d = is_all_tensor ? max_P_1d : max_P; 1173 } 1174 if (max_dim == 0) max_dim = 1; 1175 data->dim = max_dim; 1176 if (is_at_points) { 1177 CeedElemRestriction_Hip *rstr_data; 1178 CeedElemRestriction rstr_points = NULL; 1179 1180 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL)); 1181 CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points)); 1182 CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride)); 1183 CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data)); 1184 data->points.indices = (CeedInt *)rstr_data->d_offsets; 1185 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points)); 1186 } 1187 if (is_at_points) use_3d_slices = false; 1188 if (Q_1d == 0) { 1189 if (is_at_points) Q_1d = max_num_points; 1190 else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d)); 1191 } 1192 if (Q == 0) Q = Q_1d; 1193 data->Q = Q; 1194 data->Q_1d = Q_1d; 1195 1196 // Check for restriction only identity operator 1197 { 1198 bool is_identity_qf; 1199 1200 CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf)); 1201 if (is_identity_qf) { 1202 CeedEvalMode eval_mode_in, eval_mode_out; 1203 1204 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in)); 1205 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out)); 1206 CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND, 1207 "Backend does not implement restriction only identity operators"); 1208 } 1209 } 1210 1211 // Load basis source files 1212 if (!is_all_tensor) { 1213 code << "// Tensor basis source\n"; 1214 code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-templates.h>\n\n"; 1215 } 1216 if (!is_all_tensor) { 1217 code << "// Non-tensor basis source\n"; 1218 code << "#include <ceed/jit-source/hip/hip-shared-basis-nontensor-templates.h>\n\n"; 1219 } 1220 if (is_at_points) { 1221 code << "// AtPoints basis source\n"; 1222 code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-at-points-templates.h>\n\n"; 1223 } 1224 if (!is_all_tensor && !is_all_nontensor) { 1225 code << "// Tensor basis source\n"; 1226 code << "#include <ceed/jit-source/hip/hip-shared-basis-tensor-flattened-templates.h>\n\n"; 1227 } 1228 code << "// CodeGen operator source\n"; 1229 code << "#include <ceed/jit-source/hip/hip-gen-templates.h>\n\n"; 1230 1231 // Get QFunction name 1232 std::string qfunction_name(qf_data->qfunction_name); 1233 std::string operator_name; 1234 1235 operator_name = "CeedKernelHipGenOperator_" + qfunction_name; 1236 1237 // Define CEED_Q_VLA 1238 code << "\n#undef CEED_Q_VLA\n"; 1239 if (max_dim != 3 || is_at_points || use_3d_slices || !is_all_tensor) { 1240 code << "#define CEED_Q_VLA 1\n\n"; 1241 } else { 1242 code << "#define CEED_Q_VLA " << Q_1d << "\n\n"; 1243 } 1244 1245 // Add user QFunction source 1246 { 1247 const char *source_path; 1248 1249 CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path)); 1250 CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/hip/gen backend requires QFunction source code file"); 1251 1252 code << "// User QFunction source\n"; 1253 code << "#include \"" << source_path << "\"\n\n"; 1254 } 1255 1256 // Setup 1257 code << "\n// -----------------------------------------------------------------------------\n"; 1258 code << "// Operator Kernel\n"; 1259 code << "// \n"; 1260 code << "// d_[in,out]_i: CeedVector device array\n"; 1261 code << "// r_[in,out]_e_i: Element vector register\n"; 1262 code << "// r_[in,out]_q_i: Quadrature space vector register\n"; 1263 code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n"; 1264 code << "// r_[in,out]_s_i: Quadrature space slice vector register\n"; 1265 code << "// \n"; 1266 code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n"; 1267 code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n"; 1268 code << "// -----------------------------------------------------------------------------\n"; 1269 code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n"; 1270 code << "__global__ void " << operator_name 1271 << "(CeedInt num_elem, void* ctx, FieldsInt_Hip indices, Fields_Hip fields, Fields_Hip B, Fields_Hip G, CeedScalar* W, Points_Hip points) {\n"; 1272 1273 // Scratch buffers 1274 for (CeedInt i = 0; i < num_input_fields; i++) { 1275 CeedEvalMode eval_mode; 1276 1277 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 1278 if (eval_mode != CEED_EVAL_WEIGHT) { // Skip CEED_EVAL_WEIGHT 1279 code << " const CeedScalar *__restrict__ d_in_" << i << " = fields.inputs[" << i << "];\n"; 1280 } 1281 } 1282 for (CeedInt i = 0; i < num_output_fields; i++) { 1283 code << " CeedScalar *__restrict__ d_out_" << i << " = fields.outputs[" << i << "];\n"; 1284 } 1285 1286 code << " const CeedInt max_dim = " << max_dim << ";\n"; 1287 if (!is_all_tensor) { 1288 code << " const CeedInt Q = " << Q << ";\n"; 1289 } 1290 if (!is_all_nontensor) { 1291 code << " const CeedInt Q_1d = " << Q_1d << ";\n"; 1292 } 1293 if (is_at_points) { 1294 code << " const CeedInt max_num_points = " << max_num_points << ";\n"; 1295 code << " const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n"; 1296 } 1297 1298 // Shared data 1299 code << " extern __shared__ CeedScalar slice[];\n"; 1300 code << " SharedData_Hip data;\n"; 1301 code << " data.t_id_x = threadIdx.x;\n"; 1302 code << " data.t_id_y = threadIdx.y;\n"; 1303 code << " data.t_id_z = threadIdx.z;\n"; 1304 code << " data.t_id = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n"; 1305 code << " data.slice = slice + data.t_id_z*OP_T_1D" << ((!is_all_tensor || max_dim == 1) ? "" : "*OP_T_1D") << ";\n"; 1306 1307 // -- Determine input mat reuse 1308 FieldReuse_Hip input_matrix_reuse[CEED_FIELD_MAX]; 1309 1310 for (CeedInt i = 0; i < num_input_fields; i++) { 1311 input_matrix_reuse[i].index = -1; 1312 } 1313 for (CeedInt i = 0; i < num_input_fields; i++) { 1314 bool is_tensor = true; 1315 CeedEvalMode eval_mode_i; 1316 CeedBasis basis_i; 1317 1318 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode_i)); 1319 if (eval_mode_i == CEED_EVAL_WEIGHT) continue; 1320 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis_i)); 1321 CeedCallBackend(CeedBasisIsTensor(basis_i, &is_tensor)); 1322 for (CeedInt j = 0; (input_matrix_reuse[i].index == -1) && (j < i); j++) { 1323 CeedEvalMode eval_mode_j; 1324 CeedBasis basis_j; 1325 1326 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j)); 1327 if (eval_mode_j == CEED_EVAL_WEIGHT) continue; 1328 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j)); 1329 if (basis_i == basis_j) { 1330 if (is_tensor) { 1331 input_matrix_reuse[i].index = j; 1332 input_matrix_reuse[i].is_input = true; 1333 input_matrix_reuse[i].eval_mode = eval_mode_j; 1334 } else { 1335 // For non-tensor can only re-use with the same eval mode 1336 if (eval_mode_i == eval_mode_j) { 1337 input_matrix_reuse[i].index = j; 1338 input_matrix_reuse[i].is_input = true; 1339 input_matrix_reuse[i].eval_mode = eval_mode_j; 1340 } 1341 } 1342 } 1343 CeedCallBackend(CeedBasisDestroy(&basis_j)); 1344 } 1345 CeedCallBackend(CeedBasisDestroy(&basis_i)); 1346 } 1347 1348 // -- Determine output mat reuse 1349 FieldReuse_Hip output_matrix_reuse[CEED_FIELD_MAX]; 1350 1351 for (CeedInt i = 0; i < num_output_fields; i++) { 1352 output_matrix_reuse[i].index = -1; 1353 } 1354 for (CeedInt i = 0; i < num_output_fields; i++) { 1355 bool is_tensor = true; 1356 CeedEvalMode eval_mode_i; 1357 CeedBasis basis_i; 1358 1359 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode_i)); 1360 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis_i)); 1361 for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < num_input_fields); j++) { 1362 CeedEvalMode eval_mode_j; 1363 CeedBasis basis_j; 1364 1365 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[j], &eval_mode_j)); 1366 if (eval_mode_j == CEED_EVAL_WEIGHT) continue; 1367 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[j], &basis_j)); 1368 if (basis_i == basis_j) { 1369 if (is_tensor) { 1370 output_matrix_reuse[i].index = j; 1371 output_matrix_reuse[i].is_input = true; 1372 output_matrix_reuse[i].eval_mode = eval_mode_j; 1373 } else { 1374 // For non-tensor can only re-use with the same eval mode 1375 if (eval_mode_i == eval_mode_j) { 1376 output_matrix_reuse[i].index = j; 1377 output_matrix_reuse[i].is_input = true; 1378 output_matrix_reuse[i].eval_mode = eval_mode_j; 1379 } 1380 } 1381 } 1382 CeedCallBackend(CeedBasisDestroy(&basis_j)); 1383 } 1384 for (CeedInt j = 0; (output_matrix_reuse[i].index == -1) && (j < i); j++) { 1385 CeedEvalMode eval_mode_j; 1386 CeedBasis basis_j; 1387 1388 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode_j)); 1389 if (eval_mode_j == CEED_EVAL_WEIGHT) continue; 1390 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis_j)); 1391 CeedCallBackend(CeedBasisIsTensor(basis_i, &is_tensor)); 1392 if (basis_i == basis_j) { 1393 if (is_tensor) { 1394 output_matrix_reuse[i].index = j; 1395 output_matrix_reuse[i].is_input = false; 1396 output_matrix_reuse[i].eval_mode = eval_mode_j; 1397 } else { 1398 // For non-tensor can only re-use with the same eval mode 1399 if (eval_mode_i == eval_mode_j) { 1400 output_matrix_reuse[i].index = j; 1401 output_matrix_reuse[i].is_input = false; 1402 output_matrix_reuse[i].eval_mode = eval_mode_j; 1403 } 1404 } 1405 } 1406 CeedCallBackend(CeedBasisDestroy(&basis_j)); 1407 } 1408 CeedCallBackend(CeedBasisDestroy(&basis_i)); 1409 } 1410 1411 // Initialize constants, and matrices B and G 1412 code << "\n // Input field constants and basis data\n"; 1413 for (CeedInt i = 0; i < num_input_fields; i++) { 1414 CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_input_fields[i], qf_input_fields[i], input_matrix_reuse[i], max_dim, Q, 1415 Q_1d, true, is_all_tensor, is_at_points, use_3d_slices)); 1416 } 1417 code << "\n // Output field constants and basis data\n"; 1418 for (CeedInt i = 0; i < num_output_fields; i++) { 1419 CeedCallBackend(CeedOperatorBuildKernelFieldData_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], output_matrix_reuse[i], max_dim, 1420 Q, Q_1d, false, is_all_tensor, is_at_points, use_3d_slices)); 1421 } 1422 1423 // Loop over all elements 1424 code << "\n // Element loop\n"; 1425 code << " __syncthreads();\n"; 1426 code << " for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n"; 1427 1428 // -- Compute minimum buffer space needed 1429 CeedInt max_rstr_buffer_size = 1; 1430 1431 for (CeedInt i = 0; i < num_input_fields; i++) { 1432 CeedInt num_comp, elem_size; 1433 CeedElemRestriction elem_rstr; 1434 1435 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 1436 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 1437 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1438 max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? elem_size : 1)); 1439 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 1440 } 1441 for (CeedInt i = 0; i < num_output_fields; i++) { 1442 CeedInt num_comp, elem_size; 1443 CeedElemRestriction elem_rstr; 1444 1445 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 1446 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 1447 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1448 max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_all_tensor && (max_dim >= 3) ? elem_size : 1)); 1449 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 1450 } 1451 code << " // Scratch restriction buffer space\n"; 1452 code << " CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n"; 1453 1454 // -- Determine best input field processing order 1455 CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX]; 1456 1457 for (CeedInt i = 0; i < num_input_fields; i++) { 1458 field_rstr_in_buffer[i] = -1; 1459 input_field_order[i] = -1; 1460 } 1461 { 1462 bool is_ordered[CEED_FIELD_MAX]; 1463 CeedInt curr_index = 0; 1464 1465 for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false; 1466 for (CeedInt i = 0; i < num_input_fields; i++) { 1467 CeedVector vec_i; 1468 CeedElemRestriction rstr_i; 1469 1470 if (is_ordered[i]) continue; 1471 field_rstr_in_buffer[i] = i; 1472 is_ordered[i] = true; 1473 input_field_order[curr_index] = i; 1474 curr_index++; 1475 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i)); 1476 if (vec_i == CEED_VECTOR_NONE) continue; // CEED_EVAL_WEIGHT 1477 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i)); 1478 for (CeedInt j = i + 1; j < num_input_fields; j++) { 1479 CeedVector vec_j; 1480 CeedElemRestriction rstr_j; 1481 1482 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j)); 1483 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j)); 1484 if (rstr_i == rstr_j && vec_i == vec_j) { 1485 field_rstr_in_buffer[j] = i; 1486 is_ordered[j] = true; 1487 input_field_order[curr_index] = j; 1488 curr_index++; 1489 } 1490 CeedCallBackend(CeedVectorDestroy(&vec_j)); 1491 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j)); 1492 } 1493 CeedCallBackend(CeedVectorDestroy(&vec_i)); 1494 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i)); 1495 } 1496 } 1497 1498 // -- Input restriction and basis 1499 code << "\n // -- Input field restrictions and basis actions\n"; 1500 for (CeedInt i = 0; i < num_input_fields; i++) { 1501 const char *field_name; 1502 const CeedInt f = input_field_order[i]; 1503 1504 CeedCallBackend(CeedOperatorFieldGetName(op_input_fields[f], &field_name)); 1505 code << " // ---- Input field " << f << ": " << field_name << "\n"; 1506 1507 // ---- Restriction 1508 CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, f, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], max_dim, 1509 Q_1d, true, is_all_tensor, is_at_points, use_3d_slices)); 1510 1511 // ---- Basis action 1512 CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, f, op_input_fields[f], qf_input_fields[f], max_dim, Q_1d, true, is_all_tensor, 1513 is_at_points, use_3d_slices)); 1514 } 1515 1516 // -- Q function 1517 CeedCallBackend(CeedOperatorBuildKernelQFunction_Hip_gen(code, data, max_dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields, 1518 num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d, is_all_tensor, 1519 is_at_points, use_3d_slices)); 1520 1521 // -- Output basis and restriction 1522 code << "\n // -- Output field basis action and restrictions\n"; 1523 for (CeedInt i = 0; i < num_output_fields; i++) { 1524 const char *field_name; 1525 1526 CeedCallBackend(CeedOperatorFieldGetName(op_output_fields[i], &field_name)); 1527 code << " // ---- Output field " << i << ": " << field_name << "\n"; 1528 1529 // ---- Basis action 1530 CeedCallBackend(CeedOperatorBuildKernelBasis_Hip_gen(code, data, i, op_output_fields[i], qf_output_fields[i], max_dim, Q_1d, false, is_all_tensor, 1531 is_at_points, use_3d_slices)); 1532 1533 // ---- Restriction 1534 CeedCallBackend(CeedOperatorBuildKernelRestriction_Hip_gen(code, data, i, NULL, op_output_fields[i], qf_output_fields[i], max_dim, Q_1d, false, 1535 is_all_tensor, is_at_points, use_3d_slices)); 1536 } 1537 1538 // Close loop and function 1539 code << " }\n"; 1540 code << "}\n"; 1541 code << "// -----------------------------------------------------------------------------\n\n"; 1542 1543 CeedInt block_sizes[3] = {0, 0, 0}; 1544 CeedInt num_elem; 1545 1546 // Compile 1547 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 1548 CeedCallBackend(BlockGridCalculate_Hip_gen(is_all_tensor ? max_dim : 1, num_elem, data->max_P_1d, is_all_tensor ? Q_1d : Q, block_sizes)); 1549 if (is_at_points) block_sizes[2] = 1; 1550 { 1551 bool is_compile_good = false; 1552 1553 CeedCallBackend(CeedTryCompile_Hip(ceed, code.str().c_str(), &is_compile_good, &data->module, 2, "OP_T_1D", block_sizes[0], "BLOCK_SIZE", 1554 block_sizes[0] * block_sizes[1] * block_sizes[2])); 1555 if (is_compile_good) { 1556 *is_good_build = true; 1557 CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, operator_name.c_str(), &data->op)); 1558 } else { 1559 *is_good_build = false; 1560 data->use_fallback = true; 1561 } 1562 } 1563 CeedCallBackend(CeedOperatorSetSetupDone(op)); 1564 CeedCallBackend(CeedDestroy(&ceed)); 1565 CeedCallBackend(CeedQFunctionDestroy(&qf)); 1566 return CEED_ERROR_SUCCESS; 1567 } 1568 1569 //------------------------------------------------------------------------------ 1570