13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, 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 10ec3da8bcSJed Brown #include <ceed/backend.h> 11*2b730f8bSJeremy L Thompson #include <ceed/ceed.h> 129e201c85SYohann #include <ceed/jit-tools.h> 13*2b730f8bSJeremy L Thompson 147d8d0e25Snbeams #include <iostream> 157d8d0e25Snbeams #include <sstream> 16*2b730f8bSJeremy L Thompson #include <string> 17*2b730f8bSJeremy L Thompson 180d0321e0SJeremy L Thompson #include "../hip-ref/ceed-hip-ref.h" 197d8d0e25Snbeams #include "../hip-shared/ceed-hip-shared.h" 207d8d0e25Snbeams #include "../hip/ceed-hip-compile.h" 21*2b730f8bSJeremy L Thompson #include "ceed-hip-gen.h" 22b3e1519bSnbeams 23b3e1519bSnbeams //------------------------------------------------------------------------------ 24b3e1519bSnbeams // Calculate the block size used for launching the operator kernel 25b3e1519bSnbeams //------------------------------------------------------------------------------ 26*2b730f8bSJeremy 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) { 279e201c85SYohann const CeedInt thread1d = CeedIntMax(Q_1d, P_1d); 28b3e1519bSnbeams if (dim == 1) { 299e201c85SYohann CeedInt elems_per_block = 64 * thread1d > 256 ? 256 / thread1d : 64; 309e201c85SYohann elems_per_block = elems_per_block > 0 ? elems_per_block : 1; 31b3e1519bSnbeams block_sizes[0] = thread1d; 32b3e1519bSnbeams block_sizes[1] = 1; 339e201c85SYohann block_sizes[2] = elems_per_block; 34b3e1519bSnbeams } else if (dim == 2) { 359e201c85SYohann const CeedInt elems_per_block = thread1d < 4 ? 16 : 2; 36b3e1519bSnbeams block_sizes[0] = thread1d; 37b3e1519bSnbeams block_sizes[1] = thread1d; 389e201c85SYohann block_sizes[2] = elems_per_block; 39b3e1519bSnbeams } else if (dim == 3) { 409e201c85SYohann const CeedInt elems_per_block = thread1d < 6 ? 4 : (thread1d < 8 ? 2 : 1); 41b3e1519bSnbeams block_sizes[0] = thread1d; 42b3e1519bSnbeams block_sizes[1] = thread1d; 439e201c85SYohann block_sizes[2] = elems_per_block; 44b3e1519bSnbeams } 45b3e1519bSnbeams return CEED_ERROR_SUCCESS; 46b3e1519bSnbeams } 47b3e1519bSnbeams 487d8d0e25Snbeams //------------------------------------------------------------------------------ 499e201c85SYohann // Build single operator kernel 507d8d0e25Snbeams //------------------------------------------------------------------------------ 5137c3b1cfSnbeams extern "C" int CeedHipGenOperatorBuild(CeedOperator op) { 527d8d0e25Snbeams using std::ostringstream; 537d8d0e25Snbeams using std::string; 549e201c85SYohann bool is_setup_done; 55*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 569e201c85SYohann if (is_setup_done) return CEED_ERROR_SUCCESS; 577d8d0e25Snbeams Ceed ceed; 58*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 597d8d0e25Snbeams CeedOperator_Hip_gen *data; 60*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetData(op, &data)); 617d8d0e25Snbeams CeedQFunction qf; 627d8d0e25Snbeams CeedQFunction_Hip_gen *qf_data; 63*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 64*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 65e79b91d9SJeremy L Thompson CeedSize lsize; 66*2b730f8bSJeremy L Thompson CeedInt Q, P_1d = 0, Q_1d = 0, elem_size, num_input_fields, num_output_fields, num_comp, dim = 1; 67*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 689e201c85SYohann Q_1d = Q; 699e201c85SYohann CeedOperatorField *op_input_fields, *op_output_fields; 70*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 719e201c85SYohann CeedQFunctionField *qf_input_fields, *qf_output_fields; 72*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 739e201c85SYohann CeedEvalMode eval_mode; 747d8d0e25Snbeams CeedBasis basis; 757d8d0e25Snbeams CeedBasis_Hip_shared *basis_data; 767d8d0e25Snbeams CeedElemRestriction Erestrict; 777d8d0e25Snbeams CeedElemRestriction_Hip *restr_data; 787d8d0e25Snbeams 790b454692Sjeremylt // Check for restriction only identity operator 800b454692Sjeremylt bool is_identity_qf; 81*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf)); 820b454692Sjeremylt if (is_identity_qf) { 839e201c85SYohann CeedEvalMode eval_mode_in, eval_mode_out; 84*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in)); 85*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out)); 869e201c85SYohann if (eval_mode_in == CEED_EVAL_NONE && eval_mode_out == CEED_EVAL_NONE) 870b454692Sjeremylt // LCOV_EXCL_START 88*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement restriction only identity operators"); 890b454692Sjeremylt // LCOV_EXCL_STOP 900b454692Sjeremylt } 910b454692Sjeremylt 927d8d0e25Snbeams ostringstream code; 939e201c85SYohann // TODO: generalize to accept different device functions? 949e201c85SYohann { 959e201c85SYohann char *tensor_basis_kernel_path, *tensor_basis_kernel_source; 96*2b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-shared-basis-tensor-templates.h", &tensor_basis_kernel_path)); 979e201c85SYohann CeedDebug256(ceed, 2, "----- Loading Tensor Basis Kernel Source -----\n"); 98*2b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source)); 999e201c85SYohann code << tensor_basis_kernel_source; 100*2b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&tensor_basis_kernel_path)); 101*2b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&tensor_basis_kernel_source)); 1029e201c85SYohann } 1039e201c85SYohann { 1049e201c85SYohann char *hip_gen_template_path, *hip_gen_template_source; 105*2b730f8bSJeremy L Thompson CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-gen-templates.h", &hip_gen_template_path)); 1069e201c85SYohann CeedDebug256(ceed, 2, "----- Loading Hip-Gen Template Source -----\n"); 107*2b730f8bSJeremy L Thompson CeedCallBackend(CeedLoadSourceToBuffer(ceed, hip_gen_template_path, &hip_gen_template_source)); 1089e201c85SYohann code << hip_gen_template_source; 109*2b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&hip_gen_template_path)); 110*2b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&hip_gen_template_source)); 1119e201c85SYohann } 1127d8d0e25Snbeams 1139e201c85SYohann string q_function_source(qf_data->q_function_source); 1149e201c85SYohann string q_function_name(qf_data->q_function_name); 1159e201c85SYohann string operator_name; 116204bfdd7SJeremy L Thompson operator_name = "CeedKernelHipGenOperator_" + q_function_name; 1177d8d0e25Snbeams 1189e201c85SYohann // Find dim, P_1d, Q_1d 1199e201c85SYohann data->max_P_1d = 0; 1209e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 121*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 1227d8d0e25Snbeams if (basis != CEED_BASIS_COLLOCATED) { 123*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 124*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 1257d8d0e25Snbeams 1269e201c85SYohann // Collect dim, P_1d, and Q_1d 127*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 1287d8d0e25Snbeams bool isTensor; 129*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisIsTensor(basis, &isTensor)); 1307d8d0e25Snbeams if (isTensor) { 131*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 132*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 1339e201c85SYohann if (P_1d > data->max_P_1d) data->max_P_1d = P_1d; 1347d8d0e25Snbeams } else { 1357d8d0e25Snbeams // LCOV_EXCL_START 136e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 1377d8d0e25Snbeams // LCOV_EXCL_STOP 1387d8d0e25Snbeams } 1397d8d0e25Snbeams } 1407d8d0e25Snbeams } 1419e201c85SYohann // Check output bases for Q_1d, dim as well 142dc64899eSYohann // The only input basis might be CEED_BASIS_COLLOCATED 1439e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 144*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 1457d8d0e25Snbeams 1467d8d0e25Snbeams if (basis != CEED_BASIS_COLLOCATED) { 147*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 148*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1497d8d0e25Snbeams 1509e201c85SYohann // Collect Q_1d 151*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 1527d8d0e25Snbeams bool isTensor; 153*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisIsTensor(basis, &isTensor)); 1547d8d0e25Snbeams if (isTensor) { 155*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 1567d8d0e25Snbeams } else { 1577d8d0e25Snbeams // LCOV_EXCL_START 158e15f9bd0SJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 1597d8d0e25Snbeams // LCOV_EXCL_STOP 1607d8d0e25Snbeams } 1617d8d0e25Snbeams } 1627d8d0e25Snbeams } 1637d8d0e25Snbeams data->dim = dim; 1649e201c85SYohann data->Q_1d = Q_1d; 1657d8d0e25Snbeams 1669e201c85SYohann // Only use 3D collocated gradient parallelization strategy when gradient is computed 1679e201c85SYohann // TODO: put in a function? 1689e201c85SYohann bool use_collograd_parallelization = false; 1699e201c85SYohann if (dim == 3) { 1709e201c85SYohann bool was_grad_found = false; 1719e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 172*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 1739e201c85SYohann if (eval_mode == CEED_EVAL_GRAD) { 174*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 175*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 1769e201c85SYohann use_collograd_parallelization = !!basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true); 1779e201c85SYohann was_grad_found = true; 1789e201c85SYohann } 1799e201c85SYohann } 1809e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 181*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1829e201c85SYohann if (eval_mode == CEED_EVAL_GRAD) { 183*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 184*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 1859e201c85SYohann use_collograd_parallelization = !!basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true); 1869e201c85SYohann was_grad_found = true; 1879e201c85SYohann } 1889e201c85SYohann } 1897d8d0e25Snbeams } 1907d8d0e25Snbeams 1919e201c85SYohann // Define CEED_Q_VLA 1929e201c85SYohann code << "\n#undef CEED_Q_VLA\n"; 1939e201c85SYohann if (dim != 3 || use_collograd_parallelization) { 1949e201c85SYohann code << "#define CEED_Q_VLA 1\n\n"; 1959e201c85SYohann } else { 1969e201c85SYohann code << "#define CEED_Q_VLA " << Q_1d << "\n\n"; 1979e201c85SYohann } 1989e201c85SYohann 1999e201c85SYohann code << q_function_source; 2007d8d0e25Snbeams 2017d8d0e25Snbeams // Setup 2027d8d0e25Snbeams code << "\n// -----------------------------------------------------------------------------\n"; 203b3e1519bSnbeams code << "\nextern \"C\" __launch_bounds__(BLOCK_SIZE)\n"; 204*2b730f8bSJeremy L Thompson code << "__global__ void " << operator_name 205*2b730f8bSJeremy L Thompson << "(CeedInt num_elem, void* ctx, FieldsInt_Hip indices, Fields_Hip fields, Fields_Hip B, Fields_Hip G, CeedScalar* W) {\n"; 2069e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 207*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 2089e201c85SYohann if (eval_mode != CEED_EVAL_WEIGHT) { // Skip CEED_EVAL_WEIGHT 2099e201c85SYohann code << " const CeedScalar* d_u_" << i << " = fields.inputs[" << i << "];\n"; 2107d8d0e25Snbeams } 2117d8d0e25Snbeams } 2127d8d0e25Snbeams 2139e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 2149e201c85SYohann code << " CeedScalar* d_v_" << i << " = fields.outputs[" << i << "];\n"; 2157d8d0e25Snbeams } 2167d8d0e25Snbeams 2179e201c85SYohann code << " const CeedInt dim = " << dim << ";\n"; 2189e201c85SYohann code << " const CeedInt Q_1d = " << Q_1d << ";\n"; 2197d8d0e25Snbeams 2207d8d0e25Snbeams code << " HIP_DYNAMIC_SHARED( CeedScalar, slice)\n"; 2219e201c85SYohann code << " SharedData_Hip data;\n"; 2229e201c85SYohann code << " data.t_id_x = threadIdx.x;\n"; 2239e201c85SYohann code << " data.t_id_y = threadIdx.y;\n"; 2249e201c85SYohann code << " data.t_id_z = threadIdx.z;\n"; 2259e201c85SYohann code << " data.t_id = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n"; 2269e201c85SYohann code << " data.slice = slice+data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n"; 2277d8d0e25Snbeams 2287d8d0e25Snbeams code << "\n // -- Input field constants and basis data --\n"; 2297d8d0e25Snbeams // Initialize constants, and matrices B and G 2309e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 2317d8d0e25Snbeams code << " // ---- Input field " << i << " ----\n"; 2329e201c85SYohann // Get elem_size, eval_mode, num_comp 233*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict)); 234*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size)); 235*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 236*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp)); 2377d8d0e25Snbeams 2387d8d0e25Snbeams // Set field constants 2399e201c85SYohann if (eval_mode != CEED_EVAL_WEIGHT) { 240*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 2417d8d0e25Snbeams if (basis != CEED_BASIS_COLLOCATED) { 242*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 2439e201c85SYohann code << " const CeedInt P_in_" << i << " = " << P_1d << ";\n"; 2447d8d0e25Snbeams } else { 2459e201c85SYohann code << " const CeedInt P_in_" << i << " = " << Q_1d << ";\n"; 2467d8d0e25Snbeams } 2479e201c85SYohann code << " const CeedInt num_comp_in_" << i << " = " << num_comp << ";\n"; 2487d8d0e25Snbeams } 2497d8d0e25Snbeams 2507d8d0e25Snbeams // Load basis data 2519e201c85SYohann code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 2529e201c85SYohann switch (eval_mode) { 2537d8d0e25Snbeams case CEED_EVAL_NONE: 2547d8d0e25Snbeams break; 2557d8d0e25Snbeams case CEED_EVAL_INTERP: 256*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 2579e201c85SYohann data->B.inputs[i] = basis_data->d_interp_1d; 2589e201c85SYohann code << " __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n"; 2599e201c85SYohann code << " loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n"; 2607d8d0e25Snbeams break; 2617d8d0e25Snbeams case CEED_EVAL_GRAD: 262*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 2639e201c85SYohann data->B.inputs[i] = basis_data->d_interp_1d; 2649e201c85SYohann code << " __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n"; 2659e201c85SYohann code << " loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n"; 2669e201c85SYohann if (use_collograd_parallelization) { 2679e201c85SYohann data->G.inputs[i] = basis_data->d_collo_grad_1d; 2689e201c85SYohann code << " __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * Q_1d << "];\n"; 2699e201c85SYohann code << " loadMatrix<Q_1d,Q_1d>(data, G.inputs[" << i << "], s_G_in_" << i << ");\n"; 2707d8d0e25Snbeams } else { 2719e201c85SYohann bool has_collo_grad = !!basis_data->d_collo_grad_1d; 2729e201c85SYohann data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 2739e201c85SYohann code << " __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n"; 274*2b730f8bSJeremy L Thompson code << " loadMatrix<" << (has_collo_grad ? "Q_1d" : ("P_in_" + std::to_string(i))) << ",Q_1d>(data, G.inputs[" << i << "], s_G_in_" << i 275*2b730f8bSJeremy L Thompson << ");\n"; 2767d8d0e25Snbeams } 2777d8d0e25Snbeams break; 2787d8d0e25Snbeams case CEED_EVAL_WEIGHT: 2797d8d0e25Snbeams break; // No action 2807d8d0e25Snbeams case CEED_EVAL_DIV: 2817d8d0e25Snbeams break; // TODO: Not implemented 2827d8d0e25Snbeams case CEED_EVAL_CURL: 2837d8d0e25Snbeams break; // TODO: Not implemented 2847d8d0e25Snbeams } 2857d8d0e25Snbeams } 2867d8d0e25Snbeams 2877d8d0e25Snbeams code << "\n // -- Output field constants and basis data --\n"; 2889e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 2897d8d0e25Snbeams code << " // ---- Output field " << i << " ----\n"; 2909e201c85SYohann // Get elem_size, eval_mode, num_comp 291*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &Erestrict)); 292*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size)); 293*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 294*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp)); 2957d8d0e25Snbeams 2967d8d0e25Snbeams // Set field constants 297*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 2987d8d0e25Snbeams if (basis != CEED_BASIS_COLLOCATED) { 299*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 3009e201c85SYohann code << " const CeedInt P_out_" << i << " = " << P_1d << ";\n"; 3017d8d0e25Snbeams } else { 3029e201c85SYohann code << " const CeedInt P_out_" << i << " = " << Q_1d << ";\n"; 3037d8d0e25Snbeams } 3049e201c85SYohann code << " const CeedInt num_comp_out_" << i << " = " << num_comp << ";\n"; 3057d8d0e25Snbeams 3067d8d0e25Snbeams // Load basis data 3079e201c85SYohann code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 3089e201c85SYohann switch (eval_mode) { 3097d8d0e25Snbeams case CEED_EVAL_NONE: 3107d8d0e25Snbeams break; // No action 3117d8d0e25Snbeams case CEED_EVAL_INTERP: 312*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 3139e201c85SYohann data->B.outputs[i] = basis_data->d_interp_1d; 3149e201c85SYohann code << " __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n"; 3159e201c85SYohann code << " loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n"; 3167d8d0e25Snbeams break; 3177d8d0e25Snbeams case CEED_EVAL_GRAD: 318*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 3199e201c85SYohann data->B.outputs[i] = basis_data->d_interp_1d; 3209e201c85SYohann code << " __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n"; 3219e201c85SYohann code << " loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n"; 3229e201c85SYohann if (use_collograd_parallelization) { 3239e201c85SYohann data->G.outputs[i] = basis_data->d_collo_grad_1d; 3249e201c85SYohann code << " __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * Q_1d << "];\n"; 3259e201c85SYohann code << " loadMatrix<Q_1d,Q_1d>(data, G.outputs[" << i << "], s_G_out_" << i << ");\n"; 3267d8d0e25Snbeams } else { 3279e201c85SYohann bool has_collo_grad = !!basis_data->d_collo_grad_1d; 3289e201c85SYohann data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 3299e201c85SYohann code << " __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n"; 330*2b730f8bSJeremy L Thompson code << " loadMatrix<" << (has_collo_grad ? "Q_1d" : ("P_out_" + std::to_string(i))) << ",Q_1d>(data, G.outputs[" << i << "], s_G_out_" 331*2b730f8bSJeremy L Thompson << i << ");\n"; 3327d8d0e25Snbeams } 3337d8d0e25Snbeams break; 3347d8d0e25Snbeams // LCOV_EXCL_START 3357d8d0e25Snbeams case CEED_EVAL_WEIGHT: { 3367d8d0e25Snbeams Ceed ceed; 337*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 338*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 3397d8d0e25Snbeams break; // Should not occur 3407d8d0e25Snbeams } 3417d8d0e25Snbeams case CEED_EVAL_DIV: 3427d8d0e25Snbeams break; // TODO: Not implemented 3437d8d0e25Snbeams case CEED_EVAL_CURL: 3447d8d0e25Snbeams break; // TODO: Not implemented 3457d8d0e25Snbeams // LCOV_EXCL_STOP 3467d8d0e25Snbeams } 3477d8d0e25Snbeams } 3487d8d0e25Snbeams code << "\n // -- Element loop --\n"; 3497d8d0e25Snbeams code << " __syncthreads();\n"; 3509e201c85SYohann code << " for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n"; 3517d8d0e25Snbeams // Input basis apply if needed 3527d8d0e25Snbeams // Generate the correct eval mode code for each input 3537d8d0e25Snbeams code << " // -- Input field restrictions and basis actions --\n"; 3549e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 3557d8d0e25Snbeams code << " // ---- Input field " << i << " ----\n"; 3569e201c85SYohann // Get elem_size, eval_mode, num_comp 357*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict)); 358*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size)); 359*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 360*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp)); 3617d8d0e25Snbeams 3627d8d0e25Snbeams // Restriction 363*2b730f8bSJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_collograd_parallelization)) { 3649e201c85SYohann code << " CeedScalar r_u_" << i << "[num_comp_in_" << i << "*P_in_" << i << "];\n"; 3657d8d0e25Snbeams 3669e201c85SYohann bool is_strided; 367*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided)); 3689e201c85SYohann if (!is_strided) { 369*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize)); 3707d8d0e25Snbeams code << " const CeedInt lsize_in_" << i << " = " << lsize << ";\n"; 3719e201c85SYohann CeedInt comp_stride; 372*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride)); 3739e201c85SYohann code << " // CompStride: " << comp_stride << "\n"; 374*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data)); 3759e201c85SYohann data->indices.inputs[i] = restr_data->d_ind; 376*2b730f8bSJeremy L Thompson code << " readDofsOffset" << dim << "d<num_comp_in_" << i << ", " << comp_stride << ", P_in_" << i << ">(data, lsize_in_" << i 377*2b730f8bSJeremy L Thompson << ", elem, indices.inputs[" << i << "], d_u_" << i << ", r_u_" << i << ");\n"; 3787d8d0e25Snbeams } else { 3799e201c85SYohann bool has_backend_strides; 380*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides)); 3819e201c85SYohann CeedInt num_elem; 382*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem)); 3839e201c85SYohann CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 3849e201c85SYohann if (!has_backend_strides) { 385*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides)); 3867d8d0e25Snbeams } 3877d8d0e25Snbeams code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 388*2b730f8bSJeremy L Thompson code << " readDofsStrided" << dim << "d<num_comp_in_" << i << ",P_in_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2] 389*2b730f8bSJeremy L Thompson << ">(data, elem, d_u_" << i << ", r_u_" << i << ");\n"; 3907d8d0e25Snbeams } 3917d8d0e25Snbeams } 3927d8d0e25Snbeams 3937d8d0e25Snbeams // Basis action 3949e201c85SYohann code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 3959e201c85SYohann switch (eval_mode) { 3967d8d0e25Snbeams case CEED_EVAL_NONE: 3979e201c85SYohann if (!use_collograd_parallelization) { 3989e201c85SYohann code << " CeedScalar* r_t_" << i << " = r_u_" << i << ";\n"; 3997d8d0e25Snbeams } 4007d8d0e25Snbeams break; 4017d8d0e25Snbeams case CEED_EVAL_INTERP: 4029e201c85SYohann code << " CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n"; 403*2b730f8bSJeremy L Thompson code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_in_" << i << ",P_in_" << i << ",Q_1d>(data, r_u_" << i << ", s_B_in_" 404*2b730f8bSJeremy L Thompson << i << ", r_t_" << i << ");\n"; 4057d8d0e25Snbeams break; 4067d8d0e25Snbeams case CEED_EVAL_GRAD: 4079e201c85SYohann if (use_collograd_parallelization) { 4089e201c85SYohann code << " CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n"; 409*2b730f8bSJeremy L Thompson code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_in_" << i << ",P_in_" << i << ",Q_1d>(data, r_u_" << i 410*2b730f8bSJeremy L Thompson << ", s_B_in_" << i << ", r_t_" << i << ");\n"; 4117d8d0e25Snbeams } else { 4129e201c85SYohann CeedInt P_1d; 413*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 414*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 4159e201c85SYohann code << " CeedScalar r_t_" << i << "[num_comp_in_" << i << "*dim*Q_1d];\n"; 416*2b730f8bSJeremy L Thompson code << " Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_in_" << i 417*2b730f8bSJeremy L Thompson << ",P_in_" << i << ",Q_1d>(data, r_u_" << i << ", s_B_in_" << i << ", s_G_in_" << i << ", r_t_" << i << ");\n"; 4187d8d0e25Snbeams } 4197d8d0e25Snbeams break; 4207d8d0e25Snbeams case CEED_EVAL_WEIGHT: 4219e201c85SYohann code << " CeedScalar r_t_" << i << "[Q_1d];\n"; 422*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 423*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 424437930d1SJeremy L Thompson data->W = basis_data->d_q_weight_1d; 4259e201c85SYohann code << " Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<Q_1d>(data, W, r_t_" << i << ");\n"; 4267d8d0e25Snbeams break; // No action 4277d8d0e25Snbeams case CEED_EVAL_DIV: 4287d8d0e25Snbeams break; // TODO: Not implemented 4297d8d0e25Snbeams case CEED_EVAL_CURL: 4307d8d0e25Snbeams break; // TODO: Not implemented 4317d8d0e25Snbeams } 4327d8d0e25Snbeams } 4337d8d0e25Snbeams 4347d8d0e25Snbeams // Q function 4357d8d0e25Snbeams code << "\n // -- Output field setup --\n"; 4369e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 4377d8d0e25Snbeams code << "\n // ---- Output field " << i << " ----\n"; 438*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 439*2b730f8bSJeremy L Thompson if (eval_mode == CEED_EVAL_GRAD) { 4409e201c85SYohann if (use_collograd_parallelization) { 4417d8d0e25Snbeams // Accumulator for gradient slices 4429e201c85SYohann code << " CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n"; 4439e201c85SYohann code << " for (CeedInt i = 0; i < num_comp_out_" << i << "; i++) {\n"; 4449e201c85SYohann code << " for (CeedInt j = 0; j < Q_1d; ++j) {\n"; 4459e201c85SYohann code << " r_tt_" << i << "[j + i*Q_1d] = 0.0;\n"; 4467d8d0e25Snbeams code << " }\n"; 4477d8d0e25Snbeams code << " }\n"; 4487d8d0e25Snbeams } else { 4499e201c85SYohann code << " CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*dim*Q_1d];\n"; 4507d8d0e25Snbeams } 4517d8d0e25Snbeams } 452*2b730f8bSJeremy L Thompson if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) { 4539e201c85SYohann code << " CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n"; 4547d8d0e25Snbeams } 4557d8d0e25Snbeams } 4567d8d0e25Snbeams // We treat quadrature points per slice in 3d to save registers 4579e201c85SYohann if (use_collograd_parallelization) { 4589e201c85SYohann code << "\n // Note: Using planes of 3D elements\n"; 4597d8d0e25Snbeams code << "#pragma unroll\n"; 4609e201c85SYohann code << " for (CeedInt q = 0; q < Q_1d; q++) {\n"; 4617d8d0e25Snbeams code << " // -- Input fields --\n"; 4629e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 4637d8d0e25Snbeams code << " // ---- Input field " << i << " ----\n"; 4649e201c85SYohann // Get elem_size, eval_mode, num_comp 465*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 4667d8d0e25Snbeams // Basis action 4679e201c85SYohann code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 4689e201c85SYohann switch (eval_mode) { 4697d8d0e25Snbeams case CEED_EVAL_NONE: 4709e201c85SYohann code << " CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n"; 4717d8d0e25Snbeams 4729e201c85SYohann bool is_strided; 473*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &Erestrict)); 474*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided)); 4759e201c85SYohann if (!is_strided) { 476*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize)); 4777d8d0e25Snbeams code << " const CeedInt lsize_in_" << i << " = " << lsize << ";\n"; 4789e201c85SYohann CeedInt comp_stride; 479*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride)); 4809e201c85SYohann code << " // CompStride: " << comp_stride << "\n"; 481*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data)); 4829e201c85SYohann data->indices.inputs[i] = restr_data->d_ind; 483*2b730f8bSJeremy L Thompson code << " readSliceQuadsOffset" 484*2b730f8bSJeremy L Thompson << "3d<num_comp_in_" << i << ", " << comp_stride << ", Q_1d>(data, lsize_in_" << i << ", elem, q, indices.inputs[" << i << "], d_u_" 485*2b730f8bSJeremy L Thompson << i << ", r_q_" << i << ");\n"; 4867d8d0e25Snbeams } else { 487*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size)); 4889e201c85SYohann bool has_backend_strides; 489*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides)); 4909e201c85SYohann CeedInt num_elem; 491*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem)); 4929e201c85SYohann CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 4939e201c85SYohann if (!has_backend_strides) { 494*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides)); 4957d8d0e25Snbeams } 4967d8d0e25Snbeams code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 497*2b730f8bSJeremy L Thompson code << " readSliceQuadsStrided" 498*2b730f8bSJeremy L Thompson << "3d<num_comp_in_" << i 499*2b730f8bSJeremy L Thompson << ",Q_1d" 500*2b730f8bSJeremy L Thompson "," 501*2b730f8bSJeremy L Thompson << strides[0] << "," << strides[1] << "," << strides[2] << ">(data, elem, q, d_u_" << i << ", r_q_" << i << ");\n"; 5027d8d0e25Snbeams } 5037d8d0e25Snbeams break; 5047d8d0e25Snbeams case CEED_EVAL_INTERP: 5059e201c85SYohann code << " CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n"; 5069e201c85SYohann code << " for (CeedInt j = 0; j < num_comp_in_" << i << " ; ++j) {\n"; 5079e201c85SYohann code << " r_q_" << i << "[j] = r_t_" << i << "[q + j*Q_1d];\n"; 5087d8d0e25Snbeams code << " }\n"; 5097d8d0e25Snbeams break; 5107d8d0e25Snbeams case CEED_EVAL_GRAD: 5119e201c85SYohann code << " CeedScalar r_q_" << i << "[num_comp_in_" << i << "*dim];\n"; 5129e201c85SYohann code << " gradCollo3d<num_comp_in_" << i << ",Q_1d>(data, q, r_t_" << i << ", s_G_in_" << i << ", r_q_" << i << ");\n"; 5137d8d0e25Snbeams break; 5147d8d0e25Snbeams case CEED_EVAL_WEIGHT: 5159e201c85SYohann code << " CeedScalar r_q_" << i << "[1];\n"; 5169e201c85SYohann code << " r_q_" << i << "[0] = r_t_" << i << "[q];\n"; 5177d8d0e25Snbeams break; // No action 5187d8d0e25Snbeams case CEED_EVAL_DIV: 5197d8d0e25Snbeams break; // TODO: Not implemented 5207d8d0e25Snbeams case CEED_EVAL_CURL: 5217d8d0e25Snbeams break; // TODO: Not implemented 5227d8d0e25Snbeams } 5237d8d0e25Snbeams } 5247d8d0e25Snbeams code << "\n // -- Output fields --\n"; 5259e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 5267d8d0e25Snbeams code << " // ---- Output field " << i << " ----\n"; 527*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 5287d8d0e25Snbeams // Basis action 5299e201c85SYohann switch (eval_mode) { 5307d8d0e25Snbeams case CEED_EVAL_NONE: 5319e201c85SYohann code << " CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n"; 5327d8d0e25Snbeams break; // No action 5337d8d0e25Snbeams case CEED_EVAL_INTERP: 5349e201c85SYohann code << " CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n"; 5357d8d0e25Snbeams break; 5367d8d0e25Snbeams case CEED_EVAL_GRAD: 5379e201c85SYohann code << " CeedScalar r_qq_" << i << "[num_comp_out_" << i << "*dim];\n"; 5387d8d0e25Snbeams break; 5397d8d0e25Snbeams case CEED_EVAL_WEIGHT: 5407d8d0e25Snbeams break; // Should not occur 5417d8d0e25Snbeams case CEED_EVAL_DIV: 5427d8d0e25Snbeams break; // TODO: Not implemented 5437d8d0e25Snbeams case CEED_EVAL_CURL: 5447d8d0e25Snbeams break; // TODO: Not implemented 5457d8d0e25Snbeams } 5467d8d0e25Snbeams } 5477d8d0e25Snbeams } else { 5489e201c85SYohann code << "\n // Note: Using full elements\n"; 5497d8d0e25Snbeams code << " // -- Input fields --\n"; 5509e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 5517d8d0e25Snbeams code << " // ---- Input field " << i << " ----\n"; 5529e201c85SYohann code << " CeedScalar* r_q_" << i << " = r_t_" << i << ";\n"; 5537d8d0e25Snbeams } 5547d8d0e25Snbeams code << " // -- Output fields --\n"; 5559e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 5567d8d0e25Snbeams code << " // ---- Output field " << i << " ----\n"; 5579e201c85SYohann code << " CeedScalar* r_qq_" << i << " = r_tt_" << i << ";\n"; 5587d8d0e25Snbeams } 5597d8d0e25Snbeams } 5607d8d0e25Snbeams code << "\n // -- QFunction Inputs and outputs --\n"; 5619e201c85SYohann code << " CeedScalar* in[" << num_input_fields << "];\n"; 5629e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 5637d8d0e25Snbeams code << " // ---- Input field " << i << " ----\n"; 5649e201c85SYohann code << " in[" << i << "] = r_q_" << i << ";\n"; 5657d8d0e25Snbeams } 5669e201c85SYohann code << " CeedScalar* out[" << num_output_fields << "];\n"; 5679e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 5687d8d0e25Snbeams code << " // ---- Output field " << i << " ----\n"; 5699e201c85SYohann code << " out[" << i << "] = r_qq_" << i << ";\n"; 5707d8d0e25Snbeams } 5717d8d0e25Snbeams code << "\n // -- Apply QFunction --\n"; 5729e201c85SYohann code << " " << q_function_name << "(ctx, "; 5739e201c85SYohann if (dim != 3 || use_collograd_parallelization) { 5747d8d0e25Snbeams code << "1"; 5757d8d0e25Snbeams } else { 5769e201c85SYohann code << "Q_1d"; 5777d8d0e25Snbeams } 5787d8d0e25Snbeams code << ", in, out);\n"; 5799e201c85SYohann if (use_collograd_parallelization) { 5807d8d0e25Snbeams code << " // -- Output fields --\n"; 5819e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 5827d8d0e25Snbeams code << " // ---- Output field " << i << " ----\n"; 583*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 5847d8d0e25Snbeams // Basis action 5859e201c85SYohann code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 5869e201c85SYohann switch (eval_mode) { 5877d8d0e25Snbeams case CEED_EVAL_NONE: 5889e201c85SYohann code << " for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n"; 5899e201c85SYohann code << " r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n"; 5907d8d0e25Snbeams code << " }\n"; 5917d8d0e25Snbeams break; // No action 5927d8d0e25Snbeams case CEED_EVAL_INTERP: 5939e201c85SYohann code << " for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n"; 5949e201c85SYohann code << " r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n"; 5957d8d0e25Snbeams code << " }\n"; 5967d8d0e25Snbeams break; 5977d8d0e25Snbeams case CEED_EVAL_GRAD: 5989e201c85SYohann code << " gradColloTranspose3d<num_comp_out_" << i << ",Q_1d>(data, q, r_qq_" << i << ", s_G_out_" << i << ", r_tt_" << i << ");\n"; 5997d8d0e25Snbeams break; 6007d8d0e25Snbeams case CEED_EVAL_WEIGHT: 6017d8d0e25Snbeams break; // Should not occur 6027d8d0e25Snbeams case CEED_EVAL_DIV: 6037d8d0e25Snbeams break; // TODO: Not implemented 6047d8d0e25Snbeams case CEED_EVAL_CURL: 6057d8d0e25Snbeams break; // TODO: Not implemented 6067d8d0e25Snbeams } 6077d8d0e25Snbeams } 6087d8d0e25Snbeams code << " }\n"; 6097d8d0e25Snbeams } 6107d8d0e25Snbeams 6117d8d0e25Snbeams // Output basis apply if needed 6127d8d0e25Snbeams // Generate the correct eval mode code for each output 6137d8d0e25Snbeams code << "\n // -- Output field basis action and restrictions --\n"; 6149e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 6157d8d0e25Snbeams code << " // ---- Output field " << i << " ----\n"; 6169e201c85SYohann // Get elem_size, eval_mode, num_comp 617*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &Erestrict)); 618*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(Erestrict, &elem_size)); 619*2b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 620*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(Erestrict, &num_comp)); 6217d8d0e25Snbeams // Basis action 6229e201c85SYohann code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 6239e201c85SYohann switch (eval_mode) { 6247d8d0e25Snbeams case CEED_EVAL_NONE: 6259e201c85SYohann code << " CeedScalar* r_v_" << i << " = r_tt_" << i << ";\n"; 6267d8d0e25Snbeams break; // No action 6277d8d0e25Snbeams case CEED_EVAL_INTERP: 6289e201c85SYohann code << " CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n"; 629*2b730f8bSJeremy L Thompson code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i 630*2b730f8bSJeremy L Thompson << ", s_B_out_" << i << ", r_v_" << i << ");\n"; 6317d8d0e25Snbeams break; 6327d8d0e25Snbeams case CEED_EVAL_GRAD: 6339e201c85SYohann code << " CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n"; 6349e201c85SYohann if (use_collograd_parallelization) { 635*2b730f8bSJeremy L Thompson code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i 636*2b730f8bSJeremy L Thompson << ", s_B_out_" << i << ", r_v_" << i << ");\n"; 6377d8d0e25Snbeams } else { 6389e201c85SYohann CeedInt P_1d; 639*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 640*2b730f8bSJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 641*2b730f8bSJeremy L Thompson code << " GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_out_" << i 642*2b730f8bSJeremy L Thompson << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i << ", s_B_out_" << i << ", s_G_out_" << i << ", r_v_" << i << ");\n"; 6437d8d0e25Snbeams } 6447d8d0e25Snbeams break; 6457d8d0e25Snbeams // LCOV_EXCL_START 6467d8d0e25Snbeams case CEED_EVAL_WEIGHT: { 6477d8d0e25Snbeams Ceed ceed; 648*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 649*2b730f8bSJeremy L Thompson return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 6507d8d0e25Snbeams break; // Should not occur 6517d8d0e25Snbeams } 6527d8d0e25Snbeams case CEED_EVAL_DIV: 6537d8d0e25Snbeams break; // TODO: Not implemented 6547d8d0e25Snbeams case CEED_EVAL_CURL: 6557d8d0e25Snbeams break; // TODO: Not implemented 6567d8d0e25Snbeams // LCOV_EXCL_STOP 6577d8d0e25Snbeams } 6587d8d0e25Snbeams // Restriction 6599e201c85SYohann bool is_strided; 660*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionIsStrided(Erestrict, &is_strided)); 6619e201c85SYohann if (!is_strided) { 662*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(Erestrict, &lsize)); 6637d8d0e25Snbeams code << " const CeedInt lsize_out_" << i << " = " << lsize << ";\n"; 6649e201c85SYohann CeedInt comp_stride; 665*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(Erestrict, &comp_stride)); 6669e201c85SYohann code << " // CompStride: " << comp_stride << "\n"; 667*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(Erestrict, &restr_data)); 6689e201c85SYohann data->indices.outputs[i] = restr_data->d_ind; 669*2b730f8bSJeremy L Thompson code << " writeDofsOffset" << dim << "d<num_comp_out_" << i << ", " << comp_stride << ", P_out_" << i << ">(data, lsize_out_" << i 670*2b730f8bSJeremy L Thompson << ", elem, indices.outputs[" << i << "], r_v_" << i << ", d_v_" << i << ");\n"; 6717d8d0e25Snbeams } else { 6729e201c85SYohann bool has_backend_strides; 673*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(Erestrict, &has_backend_strides)); 6749e201c85SYohann CeedInt num_elem; 675*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(Erestrict, &num_elem)); 6769e201c85SYohann CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 6779e201c85SYohann if (!has_backend_strides) { 678*2b730f8bSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(Erestrict, &strides)); 6797d8d0e25Snbeams } 6807d8d0e25Snbeams code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 681*2b730f8bSJeremy L Thompson code << " writeDofsStrided" << dim << "d<num_comp_out_" << i << ",P_out_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2] 682*2b730f8bSJeremy L Thompson << ">(data, elem, r_v_" << i << ", d_v_" << i << ");\n"; 6837d8d0e25Snbeams } 6847d8d0e25Snbeams } 6857d8d0e25Snbeams 6867d8d0e25Snbeams code << " }\n"; 6877d8d0e25Snbeams code << "}\n"; 6887d8d0e25Snbeams code << "// -----------------------------------------------------------------------------\n\n"; 6897d8d0e25Snbeams 6907d8d0e25Snbeams // View kernel for debugging 69146dc0734SJeremy L Thompson CeedDebug256(ceed, 2, "Generated Operator Kernels:\n"); 6923f21f6b1SJeremy L Thompson CeedDebug(ceed, code.str().c_str()); 6937d8d0e25Snbeams 694539ec17dSJeremy L Thompson CeedInt block_sizes[3] = {0, 0, 0}; 6959e201c85SYohann CeedInt num_elem; 696*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 697*2b730f8bSJeremy L Thompson CeedCallBackend(BlockGridCalculate_Hip_gen(dim, num_elem, data->max_P_1d, Q_1d, block_sizes)); 698*2b730f8bSJeremy L Thompson CeedCallBackend(CeedCompileHip(ceed, code.str().c_str(), &data->module, 2, "T_1D", block_sizes[0], "BLOCK_SIZE", 699*2b730f8bSJeremy L Thompson block_sizes[0] * block_sizes[1] * block_sizes[2])); 700*2b730f8bSJeremy L Thompson CeedCallBackend(CeedGetKernelHip(ceed, data->module, operator_name.c_str(), &data->op)); 7017d8d0e25Snbeams 702*2b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorSetSetupDone(op)); 703e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 7047d8d0e25Snbeams } 7057d8d0e25Snbeams //------------------------------------------------------------------------------ 706