15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3241a4b83SYohann // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 5241a4b83SYohann // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 73d576824SJeremy L Thompson 8b8e71988SJeremy L Thompson #define CEED_DEBUG_COLOR 12 97df94212SJeremy L Thompson 1049aac155SJeremy L Thompson #include <ceed.h> 11ec3da8bcSJed Brown #include <ceed/backend.h> 129e201c85SYohann #include <ceed/jit-tools.h> 133d576824SJeremy L Thompson #include <cuda_runtime.h> 142b730f8bSJeremy L Thompson 15241a4b83SYohann #include <iostream> 16241a4b83SYohann #include <sstream> 17b2165e7aSSebastian Grimberg #include <string> 182b730f8bSJeremy L Thompson 190d0321e0SJeremy L Thompson #include "../cuda-ref/ceed-cuda-ref.h" 20241a4b83SYohann #include "../cuda-shared/ceed-cuda-shared.h" 2149aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h" 222b730f8bSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h" 232b730f8bSJeremy L Thompson #include "ceed-cuda-gen.h" 24241a4b83SYohann 25ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 264b3e95d5SJeremy L Thompson // Determine type of operator 274b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 284b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelData_Cuda_gen(Ceed ceed, CeedInt num_input_fields, CeedOperatorField *op_input_fields, 294b3e95d5SJeremy L Thompson CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, CeedOperatorField *op_output_fields, 304b3e95d5SJeremy L Thompson CeedQFunctionField *qf_output_fields, CeedInt *max_P_1d, CeedInt *Q_1d, CeedInt *dim, bool *is_tensor, 314b3e95d5SJeremy L Thompson bool *use_3d_slices) { 324b3e95d5SJeremy L Thompson // Find dim, P_1d, Q_1d 334b3e95d5SJeremy L Thompson *max_P_1d = 0; 344b3e95d5SJeremy L Thompson *Q_1d = 0; 354b3e95d5SJeremy L Thompson *dim = 0; 364b3e95d5SJeremy L Thompson *is_tensor = true; 374b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 384b3e95d5SJeremy L Thompson CeedBasis basis; 394b3e95d5SJeremy L Thompson 404b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 414b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 424b3e95d5SJeremy L Thompson bool is_field_tensor; 434b3e95d5SJeremy L Thompson CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0; 444b3e95d5SJeremy L Thompson 454b3e95d5SJeremy L Thompson // Collect dim, P_1d, and Q_1d 464b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 474b3e95d5SJeremy L Thompson CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 484b3e95d5SJeremy L Thompson *is_tensor = *is_tensor && is_field_tensor; 494b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 504b3e95d5SJeremy L Thompson *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 514b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 524b3e95d5SJeremy L Thompson CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 534b3e95d5SJeremy L Thompson *dim = field_dim; 544b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 554b3e95d5SJeremy L Thompson CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 564b3e95d5SJeremy L Thompson *Q_1d = field_Q_1d; 574b3e95d5SJeremy L Thompson } 58681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 594b3e95d5SJeremy L Thompson } 604b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 614b3e95d5SJeremy L Thompson CeedBasis basis; 624b3e95d5SJeremy L Thompson 634b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 644b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 654b3e95d5SJeremy L Thompson bool is_field_tensor; 664b3e95d5SJeremy L Thompson CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0; 674b3e95d5SJeremy L Thompson 684b3e95d5SJeremy L Thompson // Collect dim, P_1d, and Q_1d 694b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 704b3e95d5SJeremy L Thompson CeedCheck(is_field_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 714b3e95d5SJeremy L Thompson *is_tensor = *is_tensor && is_field_tensor; 724b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 734b3e95d5SJeremy L Thompson *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 744b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 754b3e95d5SJeremy L Thompson CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 764b3e95d5SJeremy L Thompson *dim = field_dim; 774b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 784b3e95d5SJeremy L Thompson CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 794b3e95d5SJeremy L Thompson *Q_1d = field_Q_1d; 804b3e95d5SJeremy L Thompson } 81681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 824b3e95d5SJeremy L Thompson } 834b3e95d5SJeremy L Thompson 844b3e95d5SJeremy L Thompson // Only use 3D collocated gradient parallelization strategy when gradient is computed 854b3e95d5SJeremy L Thompson *use_3d_slices = false; 864b3e95d5SJeremy L Thompson if (*dim == 3) { 874b3e95d5SJeremy L Thompson bool was_grad_found = false; 884b3e95d5SJeremy L Thompson 894b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 904b3e95d5SJeremy L Thompson CeedEvalMode eval_mode; 914b3e95d5SJeremy L Thompson 924b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 934b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_GRAD) { 944b3e95d5SJeremy L Thompson CeedBasis_Cuda_shared *basis_data; 954b3e95d5SJeremy L Thompson CeedBasis basis; 964b3e95d5SJeremy L Thompson 974b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 984b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 994b3e95d5SJeremy L Thompson *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 1004b3e95d5SJeremy L Thompson was_grad_found = true; 101681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 1024b3e95d5SJeremy L Thompson } 1034b3e95d5SJeremy L Thompson } 1044b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 1054b3e95d5SJeremy L Thompson CeedEvalMode eval_mode; 1064b3e95d5SJeremy L Thompson 1074b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1084b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_GRAD) { 1094b3e95d5SJeremy L Thompson CeedBasis_Cuda_shared *basis_data; 1104b3e95d5SJeremy L Thompson CeedBasis basis; 1114b3e95d5SJeremy L Thompson 1124b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 1134b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 1144b3e95d5SJeremy L Thompson *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 1154b3e95d5SJeremy L Thompson was_grad_found = true; 116681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 1174b3e95d5SJeremy L Thompson } 1184b3e95d5SJeremy L Thompson } 1194b3e95d5SJeremy L Thompson } 1204b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 1214b3e95d5SJeremy L Thompson } 1224b3e95d5SJeremy L Thompson 1234b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 1244b3e95d5SJeremy L Thompson // Setup fields 1254b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 1264b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedOperatorField op_field, 1274b3e95d5SJeremy L Thompson CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool use_3d_slices) { 1284b3e95d5SJeremy L Thompson std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 1294b3e95d5SJeremy L Thompson std::string P_name = "P_1d" + var_suffix, Q_name = "Q_1d"; 1304b3e95d5SJeremy L Thompson std::string option_name = (is_input ? "inputs" : "outputs"); 1314b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 1324b3e95d5SJeremy L Thompson CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 1334b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 1344b3e95d5SJeremy L Thompson CeedBasis_Cuda_shared *basis_data; 1354b3e95d5SJeremy L Thompson CeedBasis basis; 1364b3e95d5SJeremy L Thompson 1374b3e95d5SJeremy L Thompson code << " // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n"; 1384b3e95d5SJeremy L Thompson 1394b3e95d5SJeremy L Thompson // Get field data 1404b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 1414b3e95d5SJeremy L Thompson if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 1424b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1434b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 1444b3e95d5SJeremy L Thompson } 145681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 1464b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 1474b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 1484b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 1494b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 1504b3e95d5SJeremy L Thompson } 1514b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 1524b3e95d5SJeremy L Thompson 1534b3e95d5SJeremy L Thompson // Set field constants 1544b3e95d5SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) { 1554b3e95d5SJeremy L Thompson code << " const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n"; 1564b3e95d5SJeremy L Thompson code << " const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n"; 1574b3e95d5SJeremy L Thompson } 158681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 1594b3e95d5SJeremy L Thompson 1604b3e95d5SJeremy L Thompson // Load basis data 1614b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 1624b3e95d5SJeremy L Thompson switch (eval_mode) { 1634b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 1644b3e95d5SJeremy L Thompson break; 1654b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 1664b3e95d5SJeremy L Thompson if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 1674b3e95d5SJeremy L Thompson else data->B.outputs[i] = basis_data->d_interp_1d; 1684b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n"; 169*f815fac9SJeremy L Thompson code << " LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 1704b3e95d5SJeremy L Thompson break; 1714b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 1724b3e95d5SJeremy L Thompson if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 1734b3e95d5SJeremy L Thompson else data->B.outputs[i] = basis_data->d_interp_1d; 1744b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_1d * Q_1d << "];\n"; 175*f815fac9SJeremy L Thompson code << " LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 1764b3e95d5SJeremy L Thompson if (use_3d_slices) { 1774b3e95d5SJeremy L Thompson if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d; 1784b3e95d5SJeremy L Thompson else data->G.outputs[i] = basis_data->d_collo_grad_1d; 1794b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n"; 180*f815fac9SJeremy L Thompson code << " LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 1814b3e95d5SJeremy L Thompson } else { 1824b3e95d5SJeremy L Thompson bool has_collo_grad = basis_data->d_collo_grad_1d; 1834b3e95d5SJeremy L Thompson 1844b3e95d5SJeremy L Thompson if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 1854b3e95d5SJeremy L Thompson else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 1864b3e95d5SJeremy L Thompson if (has_collo_grad) { 1874b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * Q_1d << "];\n"; 188*f815fac9SJeremy L Thompson code << " LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 1894b3e95d5SJeremy L Thompson } else { 1904b3e95d5SJeremy L Thompson code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_1d * P_1d << "];\n"; 191*f815fac9SJeremy L Thompson code << " LoadMatrix<" << P_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 1924b3e95d5SJeremy L Thompson } 1934b3e95d5SJeremy L Thompson } 1944b3e95d5SJeremy L Thompson break; 1954b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 1964b3e95d5SJeremy L Thompson break; // No action 1974b3e95d5SJeremy L Thompson // LCOV_EXCL_START 1984b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 1994b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 2004b3e95d5SJeremy L Thompson break; // TODO: Not implemented 2014b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 2024b3e95d5SJeremy L Thompson } 2034b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 2044b3e95d5SJeremy L Thompson } 2054b3e95d5SJeremy L Thompson 2064b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 2074b3e95d5SJeremy L Thompson // Restriction 2084b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 2094b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim, 210e93651e5SJeremy L Thompson CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field, 211e93651e5SJeremy L Thompson CeedInt Q_1d, bool is_input, bool use_3d_slices) { 2124b3e95d5SJeremy L Thompson std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 2134b3e95d5SJeremy L Thompson std::string P_name = "P_1d" + var_suffix; 2144b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 2154b3e95d5SJeremy L Thompson CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 2164b3e95d5SJeremy L Thompson CeedSize l_size; 217*f815fac9SJeremy L Thompson CeedRestrictionType rstr_type = CEED_RESTRICTION_STANDARD; 2184b3e95d5SJeremy L Thompson CeedElemRestriction_Cuda *rstr_data; 2194b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 2204b3e95d5SJeremy L Thompson CeedBasis basis; 2214b3e95d5SJeremy L Thompson 2224b3e95d5SJeremy L Thompson // Get field data 2234b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 2244b3e95d5SJeremy L Thompson if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 225*f815fac9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type)); 2264b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 2274b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 2284b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 2294b3e95d5SJeremy L Thompson } 2304b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 2314b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 2324b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 2334b3e95d5SJeremy L Thompson } 234681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 2354b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 2364b3e95d5SJeremy L Thompson 2374b3e95d5SJeremy L Thompson // Restriction 2384b3e95d5SJeremy L Thompson if (is_input) { 2394b3e95d5SJeremy L Thompson // Input 240e93651e5SJeremy L Thompson if (field_input_buffer[i] != i) { 241e93651e5SJeremy L Thompson std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]); 242e93651e5SJeremy L Thompson 243e93651e5SJeremy L Thompson // Restriction was already done for previous input 244e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n"; 245e93651e5SJeremy L Thompson } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices)) { 246e93651e5SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE) { 247e93651e5SJeremy L Thompson // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated 2484b3e95d5SJeremy L Thompson code << " CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n"; 249e93651e5SJeremy L Thompson } else { 250e93651e5SJeremy L Thompson // Otherwise we're using the scratch space 251e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 252e93651e5SJeremy L Thompson } 253*f815fac9SJeremy L Thompson switch (rstr_type) { 254*f815fac9SJeremy L Thompson case CEED_RESTRICTION_STANDARD: { 2554b3e95d5SJeremy L Thompson CeedInt comp_stride; 2564b3e95d5SJeremy L Thompson 2574b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 2584b3e95d5SJeremy L Thompson code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 2594b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 2604b3e95d5SJeremy L Thompson code << " // CompStride: " << comp_stride << "\n"; 2614b3e95d5SJeremy L Thompson data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 262*f815fac9SJeremy L Thompson code << " ReadLVecStandard" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" 263*f815fac9SJeremy L Thompson << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n"; 264*f815fac9SJeremy L Thompson break; 265*f815fac9SJeremy L Thompson } 266*f815fac9SJeremy L Thompson case CEED_RESTRICTION_STRIDED: { 2674b3e95d5SJeremy L Thompson bool has_backend_strides; 2684b3e95d5SJeremy L Thompson CeedInt num_elem; 2694b3e95d5SJeremy L Thompson 2704b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 2714b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 2724b3e95d5SJeremy L Thompson CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 2734b3e95d5SJeremy L Thompson 2744b3e95d5SJeremy L Thompson if (!has_backend_strides) { 2754b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 2764b3e95d5SJeremy L Thompson } 2774b3e95d5SJeremy L Thompson code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 278*f815fac9SJeremy L Thompson code << " ReadLVecStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << "," 2794b3e95d5SJeremy L Thompson << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n"; 280*f815fac9SJeremy L Thompson break; 281*f815fac9SJeremy L Thompson } 282*f815fac9SJeremy L Thompson // LCOV_EXCL_START 283*f815fac9SJeremy L Thompson case CEED_RESTRICTION_ORIENTED: 284*f815fac9SJeremy L Thompson case CEED_RESTRICTION_CURL_ORIENTED: 285*f815fac9SJeremy L Thompson case CEED_RESTRICTION_POINTS: 286*f815fac9SJeremy L Thompson break; // TODO: Not implemented 287*f815fac9SJeremy L Thompson // LCOV_EXCL_STOP 2884b3e95d5SJeremy L Thompson } 2894b3e95d5SJeremy L Thompson } 2904b3e95d5SJeremy L Thompson } else { 2914b3e95d5SJeremy L Thompson // Output 292*f815fac9SJeremy L Thompson switch (rstr_type) { 293*f815fac9SJeremy L Thompson case CEED_RESTRICTION_STANDARD: { 2944b3e95d5SJeremy L Thompson CeedInt comp_stride; 2954b3e95d5SJeremy L Thompson 2964b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 2974b3e95d5SJeremy L Thompson code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 2984b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 2994b3e95d5SJeremy L Thompson code << " // CompStride: " << comp_stride << "\n"; 3004b3e95d5SJeremy L Thompson data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets; 301*f815fac9SJeremy L Thompson code << " WriteLVecStandard" << dim << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name << ">(data, l_size" 302*f815fac9SJeremy L Thompson << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n"; 303*f815fac9SJeremy L Thompson break; 304*f815fac9SJeremy L Thompson } 305*f815fac9SJeremy L Thompson case CEED_RESTRICTION_STRIDED: { 3064b3e95d5SJeremy L Thompson bool has_backend_strides; 3074b3e95d5SJeremy L Thompson CeedInt num_elem; 3084b3e95d5SJeremy L Thompson 3094b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 3104b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 3114b3e95d5SJeremy L Thompson CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 3124b3e95d5SJeremy L Thompson 3134b3e95d5SJeremy L Thompson if (!has_backend_strides) { 3144b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 3154b3e95d5SJeremy L Thompson } 3164b3e95d5SJeremy L Thompson code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 317*f815fac9SJeremy L Thompson code << " WriteLVecStrided" << dim << "d<num_comp" << var_suffix << ", " << P_name << "," << strides[0] << "," << strides[1] << "," 3184b3e95d5SJeremy L Thompson << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n"; 319*f815fac9SJeremy L Thompson break; 320*f815fac9SJeremy L Thompson } 321*f815fac9SJeremy L Thompson // LCOV_EXCL_START 322*f815fac9SJeremy L Thompson case CEED_RESTRICTION_ORIENTED: 323*f815fac9SJeremy L Thompson case CEED_RESTRICTION_CURL_ORIENTED: 324*f815fac9SJeremy L Thompson case CEED_RESTRICTION_POINTS: 325*f815fac9SJeremy L Thompson break; // TODO: Not implemented 326*f815fac9SJeremy L Thompson // LCOV_EXCL_STOP 3274b3e95d5SJeremy L Thompson } 3284b3e95d5SJeremy L Thompson } 329681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 3304b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 3314b3e95d5SJeremy L Thompson } 3324b3e95d5SJeremy L Thompson 3334b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 3344b3e95d5SJeremy L Thompson // Basis 3354b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 3364b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim, 3374b3e95d5SJeremy L Thompson CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, 3384b3e95d5SJeremy L Thompson bool use_3d_slices) { 3394b3e95d5SJeremy L Thompson std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 3404b3e95d5SJeremy L Thompson std::string P_name = "P_1d" + var_suffix, Q_name = "Q_1d"; 3414b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 3424b3e95d5SJeremy L Thompson CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 3434b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 3444b3e95d5SJeremy L Thompson CeedBasis basis; 3454b3e95d5SJeremy L Thompson 3464b3e95d5SJeremy L Thompson // Get field data 3474b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 3484b3e95d5SJeremy L Thompson if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 3494b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 3504b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 3514b3e95d5SJeremy L Thompson } 352681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 3534b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 3544b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 3554b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 3564b3e95d5SJeremy L Thompson } 3574b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 3584b3e95d5SJeremy L Thompson 3594b3e95d5SJeremy L Thompson // Basis 3604b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 3614b3e95d5SJeremy L Thompson if (is_input) { 3624b3e95d5SJeremy L Thompson switch (eval_mode) { 3634b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 3644b3e95d5SJeremy L Thompson if (!use_3d_slices) { 3654b3e95d5SJeremy L Thompson code << " CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n"; 3664b3e95d5SJeremy L Thompson } 3674b3e95d5SJeremy L Thompson break; 3684b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 3694b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 3704b3e95d5SJeremy L Thompson code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name 3714b3e95d5SJeremy L Thompson << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n"; 3724b3e95d5SJeremy L Thompson break; 3734b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 3744b3e95d5SJeremy L Thompson if (use_3d_slices) { 3754b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 3764b3e95d5SJeremy L Thompson code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", P_1d" << var_suffix << ", " << Q_name 3774b3e95d5SJeremy L Thompson << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", r_q" << var_suffix << ");\n"; 3784b3e95d5SJeremy L Thompson } else { 3794b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n"; 3804b3e95d5SJeremy L Thompson code << " Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" << var_suffix 3814b3e95d5SJeremy L Thompson << ", P_1d" << var_suffix << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix << ", r_q" 3824b3e95d5SJeremy L Thompson << var_suffix << ");\n"; 3834b3e95d5SJeremy L Thompson } 3844b3e95d5SJeremy L Thompson break; 3854b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: { 3864b3e95d5SJeremy L Thompson CeedBasis_Cuda_shared *basis_data; 3874b3e95d5SJeremy L Thompson 3884b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[" << Q_name << "];\n"; 3894b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 3904b3e95d5SJeremy L Thompson data->W = basis_data->d_q_weight_1d; 3914b3e95d5SJeremy L Thompson code << " Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n"; 3924b3e95d5SJeremy L Thompson break; 3934b3e95d5SJeremy L Thompson } 3944b3e95d5SJeremy L Thompson // LCOV_EXCL_START 3954b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 3964b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 3974b3e95d5SJeremy L Thompson break; // TODO: Not implemented 3984b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 3994b3e95d5SJeremy L Thompson } 4004b3e95d5SJeremy L Thompson } else { 4014b3e95d5SJeremy L Thompson switch (eval_mode) { 4024b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 4034b3e95d5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n"; 4044b3e95d5SJeremy L Thompson break; // No action 4054b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 406e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 4074b3e95d5SJeremy L Thompson code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name 4084b3e95d5SJeremy L Thompson << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 4094b3e95d5SJeremy L Thompson break; 4104b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 411e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 4124b3e95d5SJeremy L Thompson if (use_3d_slices) { 4134b3e95d5SJeremy L Thompson code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp" << var_suffix << ", " << P_name << ", " << Q_name 4144b3e95d5SJeremy L Thompson << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", r_e" << var_suffix << ");\n"; 4154b3e95d5SJeremy L Thompson } else { 4164b3e95d5SJeremy L Thompson code << " GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp" 4174b3e95d5SJeremy L Thompson << var_suffix << ", " << P_name << "," << Q_name << ">(data, r_q" << var_suffix << ", s_B" << var_suffix << ", s_G" << var_suffix 4184b3e95d5SJeremy L Thompson << ", r_e" << var_suffix << ");\n"; 4194b3e95d5SJeremy L Thompson } 4204b3e95d5SJeremy L Thompson break; 4214b3e95d5SJeremy L Thompson // LCOV_EXCL_START 4224b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 4234b3e95d5SJeremy L Thompson break; // Should not occur 4244b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 4254b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 4264b3e95d5SJeremy L Thompson break; // TODO: Not implemented 4274b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 4284b3e95d5SJeremy L Thompson } 4294b3e95d5SJeremy L Thompson } 430681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 4314b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 4324b3e95d5SJeremy L Thompson } 4334b3e95d5SJeremy L Thompson 4344b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 4354b3e95d5SJeremy L Thompson // QFunction 4364b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 4374b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt dim, CeedInt num_input_fields, 4384b3e95d5SJeremy L Thompson CeedOperatorField *op_input_fields, CeedQFunctionField *qf_input_fields, 4394b3e95d5SJeremy L Thompson CeedInt num_output_fields, CeedOperatorField *op_output_fields, 4404b3e95d5SJeremy L Thompson CeedQFunctionField *qf_output_fields, std::string qfunction_name, CeedInt Q_1d, 4414b3e95d5SJeremy L Thompson bool use_3d_slices) { 4424b3e95d5SJeremy L Thompson std::string Q_name = "Q_1d"; 4434b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 4444b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 4454b3e95d5SJeremy L Thompson 4464b3e95d5SJeremy L Thompson // Setup output arays 4474b3e95d5SJeremy L Thompson code << "\n // -- Output field setup\n"; 4484b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 4494b3e95d5SJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 4504b3e95d5SJeremy L Thompson 4514b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 4524b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 4534b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) { 4544b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 4554b3e95d5SJeremy L Thompson } 4564b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_GRAD) { 4574b3e95d5SJeremy L Thompson if (use_3d_slices) { 4584b3e95d5SJeremy L Thompson // Accumulator for gradient slices 4594b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 4604b3e95d5SJeremy L Thompson code << " for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n"; 4614b3e95d5SJeremy L Thompson code << " r_q" << var_suffix << "[i] = 0.0;\n"; 4624b3e95d5SJeremy L Thompson code << " }\n"; 4634b3e95d5SJeremy L Thompson } else { 4644b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << Q_name << "];\n"; 4654b3e95d5SJeremy L Thompson } 4664b3e95d5SJeremy L Thompson } 4674b3e95d5SJeremy L Thompson } 4684b3e95d5SJeremy L Thompson 4694b3e95d5SJeremy L Thompson // We treat quadrature points per slice in 3d to save registers 4704b3e95d5SJeremy L Thompson if (use_3d_slices) { 4714b3e95d5SJeremy L Thompson code << "\n // Note: Using planes of 3D elements\n"; 4724b3e95d5SJeremy L Thompson code << " #pragma unroll\n"; 4734b3e95d5SJeremy L Thompson code << " for (CeedInt q = 0; q < " << Q_name << "; q++) {\n"; 4744b3e95d5SJeremy L Thompson code << " // -- Input fields\n"; 4754b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 4764b3e95d5SJeremy L Thompson std::string var_suffix = "_in_" + std::to_string(i); 4774b3e95d5SJeremy L Thompson 4784b3e95d5SJeremy L Thompson code << " // ---- Input field " << i << "\n"; 4794b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 4804b3e95d5SJeremy L Thompson // Basis action 4814b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 4824b3e95d5SJeremy L Thompson switch (eval_mode) { 4834b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 4844b3e95d5SJeremy L Thompson bool is_strided; 4854b3e95d5SJeremy L Thompson 4864b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 4874b3e95d5SJeremy L Thompson 4884b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 4894b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 4904b3e95d5SJeremy L Thompson if (is_strided) { 4914b3e95d5SJeremy L Thompson bool has_backend_strides; 4924b3e95d5SJeremy L Thompson CeedInt num_elem, elem_size; 4934b3e95d5SJeremy L Thompson 4944b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 4954b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 4964b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 4974b3e95d5SJeremy L Thompson CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 4984b3e95d5SJeremy L Thompson 4994b3e95d5SJeremy L Thompson if (!has_backend_strides) { 5004b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 5014b3e95d5SJeremy L Thompson } 5024b3e95d5SJeremy L Thompson code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 503*f815fac9SJeremy L Thompson code << " ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << "," << strides[0] << "," << strides[1] << "," 5044b3e95d5SJeremy L Thompson << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n"; 5054b3e95d5SJeremy L Thompson } else { 5064b3e95d5SJeremy L Thompson CeedSize l_size = 0; 5074b3e95d5SJeremy L Thompson CeedInt comp_stride; 5084b3e95d5SJeremy L Thompson CeedElemRestriction_Cuda *rstr_data; 5094b3e95d5SJeremy L Thompson 5104b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 5114b3e95d5SJeremy L Thompson code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 5124b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 5134b3e95d5SJeremy L Thompson code << " // CompStride: " << comp_stride << "\n"; 5144b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 5154b3e95d5SJeremy L Thompson data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 516*f815fac9SJeremy L Thompson code << " ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix 5174b3e95d5SJeremy L Thompson << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n"; 5184b3e95d5SJeremy L Thompson } 519681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 5204b3e95d5SJeremy L Thompson break; 5214b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 5224b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 5234b3e95d5SJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n"; 5244b3e95d5SJeremy L Thompson code << " r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n"; 5254b3e95d5SJeremy L Thompson code << " }\n"; 5264b3e95d5SJeremy L Thompson break; 5274b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 5284b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 529*f815fac9SJeremy L Thompson code << " GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix 530*f815fac9SJeremy L Thompson << ", r_s" << var_suffix << ");\n"; 5314b3e95d5SJeremy L Thompson break; 5324b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 5334b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[1];\n"; 5344b3e95d5SJeremy L Thompson code << " r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n"; 5354b3e95d5SJeremy L Thompson break; // No action 5364b3e95d5SJeremy L Thompson // LCOV_EXCL_START 5374b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 5384b3e95d5SJeremy L Thompson break; // TODO: Not implemented 5394b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 5404b3e95d5SJeremy L Thompson break; // TODO: Not implemented 5414b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 5424b3e95d5SJeremy L Thompson } 5434b3e95d5SJeremy L Thompson } 5444b3e95d5SJeremy L Thompson code << "\n // -- Output fields\n"; 5454b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 5464b3e95d5SJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 5474b3e95d5SJeremy L Thompson 5484b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 5494b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 5504b3e95d5SJeremy L Thompson // Basis action 5514b3e95d5SJeremy L Thompson switch (eval_mode) { 5524b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 5534b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 5544b3e95d5SJeremy L Thompson break; // No action 5554b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 5564b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 5574b3e95d5SJeremy L Thompson break; 5584b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 5594b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 5604b3e95d5SJeremy L Thompson break; 5614b3e95d5SJeremy L Thompson // LCOV_EXCL_START 5624b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 5634b3e95d5SJeremy L Thompson break; // Should not occur 5644b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 5654b3e95d5SJeremy L Thompson break; // TODO: Not implemented 5664b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 5674b3e95d5SJeremy L Thompson break; // TODO: Not implemented 5684b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 5694b3e95d5SJeremy L Thompson } 5704b3e95d5SJeremy L Thompson } 5714b3e95d5SJeremy L Thompson } else { 5724b3e95d5SJeremy L Thompson code << "\n // Note: Using full elements\n"; 5734b3e95d5SJeremy L Thompson code << " {\n"; 5744b3e95d5SJeremy L Thompson code << " // -- Input fields\n"; 5754b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 5764b3e95d5SJeremy L Thompson code << " // ---- Input field " << i << "\n"; 5774b3e95d5SJeremy L Thompson code << " CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n"; 5784b3e95d5SJeremy L Thompson } 5794b3e95d5SJeremy L Thompson code << " // -- Output fields\n"; 5804b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 5814b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 5824b3e95d5SJeremy L Thompson code << " CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n"; 5834b3e95d5SJeremy L Thompson } 5844b3e95d5SJeremy L Thompson } 5854b3e95d5SJeremy L Thompson 5864b3e95d5SJeremy L Thompson // Input and output buffers 5874b3e95d5SJeremy L Thompson code << "\n // -- QFunction inputs and outputs\n"; 5884b3e95d5SJeremy L Thompson code << " // ---- Inputs\n"; 5894b3e95d5SJeremy L Thompson code << " CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n"; 5904b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 5914b3e95d5SJeremy L Thompson code << " // ------ Input field " << i << "\n"; 5924b3e95d5SJeremy L Thompson code << " inputs[" << i << "] = r_s_in_" << i << ";\n"; 5934b3e95d5SJeremy L Thompson } 5944b3e95d5SJeremy L Thompson code << " // ---- Outputs\n"; 5954b3e95d5SJeremy L Thompson code << " CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n"; 5964b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 5974b3e95d5SJeremy L Thompson code << " // ------ Output field " << i << "\n"; 5984b3e95d5SJeremy L Thompson code << " outputs[" << i << "] = r_s_out_" << i << ";\n"; 5994b3e95d5SJeremy L Thompson } 6004b3e95d5SJeremy L Thompson 6014b3e95d5SJeremy L Thompson // Apply QFunction 6024b3e95d5SJeremy L Thompson code << "\n // -- Apply QFunction\n"; 6034b3e95d5SJeremy L Thompson code << " " << qfunction_name << "(ctx, "; 6044b3e95d5SJeremy L Thompson if (dim != 3 || use_3d_slices) { 6054b3e95d5SJeremy L Thompson code << "1"; 6064b3e95d5SJeremy L Thompson } else { 6074b3e95d5SJeremy L Thompson code << "Q_1d"; 6084b3e95d5SJeremy L Thompson } 6094b3e95d5SJeremy L Thompson code << ", inputs, outputs);\n"; 6104b3e95d5SJeremy L Thompson 6114b3e95d5SJeremy L Thompson // Copy or apply transpose grad, if needed 6124b3e95d5SJeremy L Thompson if (use_3d_slices) { 6134b3e95d5SJeremy L Thompson code << " // -- Output fields\n"; 6144b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 6154b3e95d5SJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 6164b3e95d5SJeremy L Thompson std::string P_name = "P_1d" + var_suffix; 6174b3e95d5SJeremy L Thompson 6184b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 6194b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 6204b3e95d5SJeremy L Thompson // Basis action 6214b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 6224b3e95d5SJeremy L Thompson switch (eval_mode) { 6234b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 6244b3e95d5SJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 6254b3e95d5SJeremy L Thompson code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 6264b3e95d5SJeremy L Thompson code << " }\n"; 6274b3e95d5SJeremy L Thompson break; // No action 6284b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 6294b3e95d5SJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 6304b3e95d5SJeremy L Thompson code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 6314b3e95d5SJeremy L Thompson code << " }\n"; 6324b3e95d5SJeremy L Thompson break; 6334b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 634*f815fac9SJeremy L Thompson code << " GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G" 635*f815fac9SJeremy L Thompson << var_suffix << ", r_q" << var_suffix << ");\n"; 6364b3e95d5SJeremy L Thompson break; 6374b3e95d5SJeremy L Thompson // LCOV_EXCL_START 6384b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 6394b3e95d5SJeremy L Thompson break; // Should not occur 6404b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 6414b3e95d5SJeremy L Thompson break; // TODO: Not implemented 6424b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 6434b3e95d5SJeremy L Thompson break; // TODO: Not implemented 6444b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 6454b3e95d5SJeremy L Thompson } 6464b3e95d5SJeremy L Thompson } 6474b3e95d5SJeremy L Thompson } 6484b3e95d5SJeremy L Thompson code << " }\n"; 6494b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 6504b3e95d5SJeremy L Thompson } 6514b3e95d5SJeremy L Thompson 6524b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 653b2165e7aSSebastian Grimberg // Build single operator kernel 654ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 655eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op) { 6564b3e95d5SJeremy L Thompson bool is_tensor = true, use_3d_slices = false; 657241a4b83SYohann Ceed ceed; 6584b3e95d5SJeremy L Thompson CeedInt Q_1d, num_input_fields, num_output_fields, dim = 1; 659ca735530SJeremy L Thompson CeedQFunctionField *qf_input_fields, *qf_output_fields; 660ca735530SJeremy L Thompson CeedQFunction_Cuda_gen *qf_data; 661ca735530SJeremy L Thompson CeedQFunction qf; 662ca735530SJeremy L Thompson CeedOperatorField *op_input_fields, *op_output_fields; 663ca735530SJeremy L Thompson CeedOperator_Cuda_gen *data; 6644b3e95d5SJeremy L Thompson std::ostringstream code; 6654b3e95d5SJeremy L Thompson 6664b3e95d5SJeremy L Thompson { 6674b3e95d5SJeremy L Thompson bool is_setup_done; 668ca735530SJeremy L Thompson 669ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 670ca735530SJeremy L Thompson if (is_setup_done) return CEED_ERROR_SUCCESS; 6714b3e95d5SJeremy L Thompson } 672ca735530SJeremy L Thompson 673ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 674ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorGetData(op, &data)); 675ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 676ca735530SJeremy L Thompson CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 677ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 678ca735530SJeremy L Thompson CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 679241a4b83SYohann 6804b3e95d5SJeremy L Thompson // Get operator data 6814b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelData_Cuda_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields, 6824b3e95d5SJeremy L Thompson qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices)); 6834b3e95d5SJeremy L Thompson if (dim == 0) dim = 1; 6844b3e95d5SJeremy L Thompson data->dim = dim; 6854b3e95d5SJeremy L Thompson if (Q_1d == 0) { 6864b3e95d5SJeremy L Thompson CeedInt Q; 6874b3e95d5SJeremy L Thompson 6884b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 6894b3e95d5SJeremy L Thompson Q_1d = Q; 6904b3e95d5SJeremy L Thompson } 6914b3e95d5SJeremy L Thompson data->Q_1d = Q_1d; 6924b3e95d5SJeremy L Thompson 6930b454692Sjeremylt // Check for restriction only identity operator 6944b3e95d5SJeremy L Thompson { 6954b3e95d5SJeremy L Thompson bool is_identity_qf; 6964b3e95d5SJeremy L Thompson 6972b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf)); 6980b454692Sjeremylt if (is_identity_qf) { 6999e201c85SYohann CeedEvalMode eval_mode_in, eval_mode_out; 700ca735530SJeremy L Thompson 7012b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in)); 7022b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out)); 7036574a04fSJeremy L Thompson CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND, 7046574a04fSJeremy L Thompson "Backend does not implement restriction only identity operators"); 7050b454692Sjeremylt } 7064b3e95d5SJeremy L Thompson } 7070b454692Sjeremylt 708f1a13f77SYohann Dudouit // Add atomicAdd function for old NVidia architectures 7094b3e95d5SJeremy L Thompson { 7104b3e95d5SJeremy L Thompson Ceed_Cuda *ceed_data; 7114b3e95d5SJeremy L Thompson struct cudaDeviceProp prop; 7124b3e95d5SJeremy L Thompson 7132b730f8bSJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_data)); 7142b730f8bSJeremy L Thompson CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id)); 71580a9ef05SNatalie Beams if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) { 7169c25dd66SJeremy L Thompson code << "// AtomicAdd fallback source\n"; 7179c25dd66SJeremy L Thompson code << "#include <ceed/jit-source/cuda/cuda-atomic-add-fallback.h>\n\n"; 718f1a13f77SYohann Dudouit } 7194b3e95d5SJeremy L Thompson } 720f1a13f77SYohann Dudouit 7219e201c85SYohann // Load basis source files 7224b3e95d5SJeremy L Thompson // TODO: Add non-tensor, AtPoints 7239c25dd66SJeremy L Thompson code << "// Tensor basis source\n"; 7249c25dd66SJeremy L Thompson code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h>\n\n"; 7259c25dd66SJeremy L Thompson code << "// CodeGen operator source\n"; 7269c25dd66SJeremy L Thompson code << "#include <ceed/jit-source/cuda/cuda-gen-templates.h>\n\n"; 727241a4b83SYohann 7284b3e95d5SJeremy L Thompson // Get QFunction name 7294b3e95d5SJeremy L Thompson std::string qfunction_name(qf_data->qfunction_name); 7304b3e95d5SJeremy L Thompson std::string operator_name; 7314b3e95d5SJeremy L Thompson 73209095acaSJeremy L Thompson operator_name = "CeedKernelCudaGenOperator_" + qfunction_name; 7334d537eeaSYohann 7349e201c85SYohann // Define CEED_Q_VLA 7359e201c85SYohann code << "\n#undef CEED_Q_VLA\n"; 7364b3e95d5SJeremy L Thompson if (dim != 3 || use_3d_slices) { 7379e201c85SYohann code << "#define CEED_Q_VLA 1\n\n"; 7389e201c85SYohann } else { 7399e201c85SYohann code << "#define CEED_Q_VLA " << Q_1d << "\n\n"; 7409e201c85SYohann } 7419e201c85SYohann 7424b3e95d5SJeremy L Thompson // Add user QFunction source 7434b3e95d5SJeremy L Thompson { 7449c25dd66SJeremy L Thompson const char *source_path; 7454b3e95d5SJeremy L Thompson 7469c25dd66SJeremy L Thompson CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path)); 7479c25dd66SJeremy L Thompson CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/cuda/gen backend requires QFunction source code file"); 7489c25dd66SJeremy L Thompson 7499c25dd66SJeremy L Thompson code << "// User QFunction source\n"; 7509c25dd66SJeremy L Thompson code << "#include \"" << source_path << "\"\n\n"; 7514b3e95d5SJeremy L Thompson } 752241a4b83SYohann 753241a4b83SYohann // Setup 754818e0025SJeremy L Thompson code << "\n// -----------------------------------------------------------------------------\n"; 7554b3e95d5SJeremy L Thompson code << "// Operator Kernel\n"; 7564b3e95d5SJeremy L Thompson code << "// \n"; 7574b3e95d5SJeremy L Thompson code << "// d_[in,out]_i: CeedVector device array\n"; 7584b3e95d5SJeremy L Thompson code << "// r_[in,out]_e_i: Element vector register\n"; 7594b3e95d5SJeremy L Thompson code << "// r_[in,out]_q_i: Quadrature space vector register\n"; 7604b3e95d5SJeremy L Thompson code << "// r_[in,out]_s_i: Quadrature space slice vector register\n"; 7614b3e95d5SJeremy L Thompson code << "// \n"; 7624b3e95d5SJeremy L Thompson code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n"; 7634b3e95d5SJeremy L Thompson code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n"; 7644b3e95d5SJeremy L Thompson code << "// -----------------------------------------------------------------------------\n"; 7654b3e95d5SJeremy L Thompson code << "extern \"C\" __global__ void " << operator_name 7662b730f8bSJeremy L Thompson << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W) {\n"; 7674b3e95d5SJeremy L Thompson 7684b3e95d5SJeremy L Thompson // Scratch buffers 7699e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 7704b3e95d5SJeremy L Thompson CeedEvalMode eval_mode; 7714b3e95d5SJeremy L Thompson 7722b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 7739e201c85SYohann if (eval_mode != CEED_EVAL_WEIGHT) { // Skip CEED_EVAL_WEIGHT 7744b3e95d5SJeremy L Thompson code << " const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n"; 775241a4b83SYohann } 776241a4b83SYohann } 7779e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 7784b3e95d5SJeremy L Thompson code << " CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n"; 779241a4b83SYohann } 7801da99368SJeremy L Thompson 7819e201c85SYohann code << " const CeedInt dim = " << dim << ";\n"; 7829e201c85SYohann code << " const CeedInt Q_1d = " << Q_1d << ";\n"; 7831da99368SJeremy L Thompson 7844b3e95d5SJeremy L Thompson // Shared data 785241a4b83SYohann code << " extern __shared__ CeedScalar slice[];\n"; 7869e201c85SYohann code << " SharedData_Cuda data;\n"; 7879e201c85SYohann code << " data.t_id_x = threadIdx.x;\n"; 7889e201c85SYohann code << " data.t_id_y = threadIdx.y;\n"; 7899e201c85SYohann code << " data.t_id_z = threadIdx.z;\n"; 7909e201c85SYohann code << " data.t_id = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n"; 7919e201c85SYohann code << " data.slice = slice + data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n"; 792920dcdc4Sjeremylt 793ac421f39SYohann // Initialize constants, and matrices B and G 7944b3e95d5SJeremy L Thompson code << "\n // Input field constants and basis data\n"; 7959e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 7965a5594ffSJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_input_fields[i], qf_input_fields[i], Q_1d, true, use_3d_slices)); 797920dcdc4Sjeremylt } 7984b3e95d5SJeremy L Thompson code << "\n // Output field constants and basis data\n"; 7999e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 8005a5594ffSJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices)); 8014b3e95d5SJeremy L Thompson } 802920dcdc4Sjeremylt 8034b3e95d5SJeremy L Thompson // Loop over all elements 8044b3e95d5SJeremy L Thompson code << "\n // Element loop\n"; 805ac421f39SYohann code << " __syncthreads();\n"; 8069e201c85SYohann code << " for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n"; 8074b3e95d5SJeremy L Thompson 808e93651e5SJeremy L Thompson // -- Compute minimum buffer space needed 809e93651e5SJeremy L Thompson CeedInt max_rstr_buffer_size = 0; 810e93651e5SJeremy L Thompson 8119e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 812e93651e5SJeremy L Thompson CeedInt num_comp, elem_size; 813e93651e5SJeremy L Thompson CeedElemRestriction elem_rstr; 814e93651e5SJeremy L Thompson 815e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 816e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 817e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 818e93651e5SJeremy L Thompson max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size); 819681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 820e93651e5SJeremy L Thompson } 821e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 822e93651e5SJeremy L Thompson CeedInt num_comp, elem_size; 823e93651e5SJeremy L Thompson CeedElemRestriction elem_rstr; 824e93651e5SJeremy L Thompson 825e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 826e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 827e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 828e93651e5SJeremy L Thompson max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * elem_size); 829681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 830e93651e5SJeremy L Thompson } 831e93651e5SJeremy L Thompson code << " // Scratch restriction buffer space\n"; 832e93651e5SJeremy L Thompson code << " CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n"; 833e93651e5SJeremy L Thompson 834e93651e5SJeremy L Thompson // -- Determine best input field processing order 835e93651e5SJeremy L Thompson CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX]; 836e93651e5SJeremy L Thompson 837e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 838e93651e5SJeremy L Thompson field_rstr_in_buffer[i] = -1; 839e93651e5SJeremy L Thompson input_field_order[i] = -1; 840e93651e5SJeremy L Thompson } 841e93651e5SJeremy L Thompson { 842e93651e5SJeremy L Thompson bool is_ordered[CEED_FIELD_MAX]; 843e93651e5SJeremy L Thompson CeedInt curr_index = 0; 844e93651e5SJeremy L Thompson 845e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false; 846e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 847e93651e5SJeremy L Thompson CeedVector vec_i; 848e93651e5SJeremy L Thompson CeedElemRestriction rstr_i; 849e93651e5SJeremy L Thompson 850e93651e5SJeremy L Thompson if (is_ordered[i]) continue; 851e93651e5SJeremy L Thompson field_rstr_in_buffer[i] = i; 852e93651e5SJeremy L Thompson is_ordered[i] = true; 853e93651e5SJeremy L Thompson input_field_order[curr_index] = i; 854e93651e5SJeremy L Thompson curr_index++; 855034f99fdSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i)); 856e93651e5SJeremy L Thompson if (vec_i == CEED_VECTOR_NONE) continue; // CEED_EVAL_WEIGHT 857e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i)); 858e93651e5SJeremy L Thompson for (CeedInt j = i + 1; j < num_input_fields; j++) { 859e93651e5SJeremy L Thompson CeedVector vec_j; 860e93651e5SJeremy L Thompson CeedElemRestriction rstr_j; 861e93651e5SJeremy L Thompson 862e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j)); 863e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j)); 864e93651e5SJeremy L Thompson if (rstr_i == rstr_j && vec_i == vec_j) { 865e93651e5SJeremy L Thompson field_rstr_in_buffer[j] = i; 866e93651e5SJeremy L Thompson is_ordered[j] = true; 867e93651e5SJeremy L Thompson input_field_order[curr_index] = j; 868e93651e5SJeremy L Thompson curr_index++; 869e93651e5SJeremy L Thompson } 870681d0ea7SJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec_j)); 871681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j)); 872e93651e5SJeremy L Thompson } 873681d0ea7SJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec_i)); 874681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i)); 875e93651e5SJeremy L Thompson } 876e93651e5SJeremy L Thompson } 877e93651e5SJeremy L Thompson 878e93651e5SJeremy L Thompson // -- Input restriction and basis 879e93651e5SJeremy L Thompson code << "\n // -- Input field restrictions and basis actions\n"; 880e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 881e93651e5SJeremy L Thompson CeedInt f = input_field_order[i]; 882e93651e5SJeremy L Thompson 883e93651e5SJeremy L Thompson code << " // ---- Input field " << f << "\n"; 884920dcdc4Sjeremylt 8854b3e95d5SJeremy L Thompson // ---- Restriction 886e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], 887e93651e5SJeremy L Thompson Q_1d, true, use_3d_slices)); 8889a2291e3SJeremy L Thompson 8894b3e95d5SJeremy L Thompson // ---- Basis action 890e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, use_3d_slices)); 891920dcdc4Sjeremylt } 892920dcdc4Sjeremylt 8934b3e95d5SJeremy L Thompson // -- Q function 8944b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, dim, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, 8954b3e95d5SJeremy L Thompson op_output_fields, qf_output_fields, qfunction_name, Q_1d, use_3d_slices)); 896ca735530SJeremy L Thompson 8974b3e95d5SJeremy L Thompson // -- Output basis and restriction 8984b3e95d5SJeremy L Thompson code << "\n // -- Output field basis action and restrictions\n"; 8999e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 9004b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 901ca735530SJeremy L Thompson 9024b3e95d5SJeremy L Thompson // ---- Basis action 9034b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices)); 9049a2291e3SJeremy L Thompson 9054b3e95d5SJeremy L Thompson // ---- Restriction 9064b3e95d5SJeremy L Thompson CeedCallBackend( 907e93651e5SJeremy L Thompson CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false, use_3d_slices)); 908ac421f39SYohann } 909241a4b83SYohann 9104b3e95d5SJeremy L Thompson // Close loop and function 911241a4b83SYohann code << " }\n"; 912818e0025SJeremy L Thompson code << "}\n"; 913818e0025SJeremy L Thompson code << "// -----------------------------------------------------------------------------\n\n"; 914241a4b83SYohann 915ab213215SJeremy L Thompson // View kernel for debugging 916eb7e6cafSJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, code.str().c_str(), &data->module, 1, "T_1D", CeedIntMax(Q_1d, data->max_P_1d))); 917eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op)); 9182b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorSetSetupDone(op)); 9199bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 920c11e12f4SJeremy L Thompson CeedCallBackend(CeedQFunctionDestroy(&qf)); 921e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 922241a4b83SYohann } 9232a86cc9dSSebastian Grimberg 924ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 925