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; 37dc007f05SJeremy L Thompson 384b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 394b3e95d5SJeremy L Thompson CeedBasis basis; 404b3e95d5SJeremy L Thompson 414b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 424b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 434b3e95d5SJeremy L Thompson bool is_field_tensor; 444b3e95d5SJeremy L Thompson CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0; 454b3e95d5SJeremy L Thompson 464b3e95d5SJeremy L Thompson // Collect dim, P_1d, and Q_1d 474b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 484b3e95d5SJeremy L Thompson *is_tensor = *is_tensor && is_field_tensor; 49dc007f05SJeremy L Thompson if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 50dc007f05SJeremy L Thompson else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d)); 514b3e95d5SJeremy L Thompson *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 524b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 534b3e95d5SJeremy L Thompson CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 544b3e95d5SJeremy L Thompson *dim = field_dim; 55dc007f05SJeremy L Thompson if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 56dc007f05SJeremy L Thompson else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d)); 574b3e95d5SJeremy L Thompson CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 584b3e95d5SJeremy L Thompson *Q_1d = field_Q_1d; 594b3e95d5SJeremy L Thompson } 60681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 614b3e95d5SJeremy L Thompson } 624b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 634b3e95d5SJeremy L Thompson CeedBasis basis; 644b3e95d5SJeremy L Thompson 654b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 664b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 674b3e95d5SJeremy L Thompson bool is_field_tensor; 684b3e95d5SJeremy L Thompson CeedInt field_P_1d = 0, field_Q_1d = 0, field_dim = 0; 694b3e95d5SJeremy L Thompson 704b3e95d5SJeremy L Thompson // Collect dim, P_1d, and Q_1d 714b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisIsTensor(basis, &is_field_tensor)); 724b3e95d5SJeremy L Thompson *is_tensor = *is_tensor && is_field_tensor; 73dc007f05SJeremy L Thompson if (is_field_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &field_P_1d)); 74dc007f05SJeremy L Thompson else CeedCallBackend(CeedBasisGetNumNodes(basis, &field_P_1d)); 754b3e95d5SJeremy L Thompson *max_P_1d = CeedIntMax(*max_P_1d, field_P_1d); 764b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetDimension(basis, &field_dim)); 774b3e95d5SJeremy L Thompson CeedCheck(*dim == 0 || field_dim == *dim, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 784b3e95d5SJeremy L Thompson *dim = field_dim; 79dc007f05SJeremy L Thompson if (is_field_tensor) CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &field_Q_1d)); 80dc007f05SJeremy L Thompson else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis, &field_Q_1d)); 814b3e95d5SJeremy L Thompson CeedCheck(*Q_1d == 0 || field_Q_1d == *Q_1d, ceed, CEED_ERROR_BACKEND, "Quadrature spaces must be compatible"); 824b3e95d5SJeremy L Thompson *Q_1d = field_Q_1d; 834b3e95d5SJeremy L Thompson } 84681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 854b3e95d5SJeremy L Thompson } 864b3e95d5SJeremy L Thompson 874b3e95d5SJeremy L Thompson // Only use 3D collocated gradient parallelization strategy when gradient is computed 884b3e95d5SJeremy L Thompson *use_3d_slices = false; 894b3e95d5SJeremy L Thompson if (*dim == 3) { 904b3e95d5SJeremy L Thompson bool was_grad_found = false; 914b3e95d5SJeremy L Thompson 924b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 934b3e95d5SJeremy L Thompson CeedEvalMode eval_mode; 944b3e95d5SJeremy L Thompson 954b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 964b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_GRAD) { 974b3e95d5SJeremy L Thompson CeedBasis_Cuda_shared *basis_data; 984b3e95d5SJeremy L Thompson CeedBasis basis; 994b3e95d5SJeremy L Thompson 1004b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 1014b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 1024b3e95d5SJeremy L Thompson *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 1034b3e95d5SJeremy L Thompson was_grad_found = true; 104681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 1054b3e95d5SJeremy L Thompson } 1064b3e95d5SJeremy L Thompson } 1074b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 1084b3e95d5SJeremy L Thompson CeedEvalMode eval_mode; 1094b3e95d5SJeremy L Thompson 1104b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1114b3e95d5SJeremy L Thompson if (eval_mode == CEED_EVAL_GRAD) { 1124b3e95d5SJeremy L Thompson CeedBasis_Cuda_shared *basis_data; 1134b3e95d5SJeremy L Thompson CeedBasis basis; 1144b3e95d5SJeremy L Thompson 1154b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 1164b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 1174b3e95d5SJeremy L Thompson *use_3d_slices = basis_data->d_collo_grad_1d && (was_grad_found ? *use_3d_slices : true); 1184b3e95d5SJeremy L Thompson was_grad_found = true; 119681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 1204b3e95d5SJeremy L Thompson } 1214b3e95d5SJeremy L Thompson } 1224b3e95d5SJeremy L Thompson } 1234b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 1244b3e95d5SJeremy L Thompson } 1254b3e95d5SJeremy L Thompson 1264b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 1274b3e95d5SJeremy L Thompson // Setup fields 1284b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 1294b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelFieldData_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedOperatorField op_field, 130dc007f05SJeremy L Thompson CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_tensor, bool is_at_points, 1318b97b69aSJeremy L Thompson bool use_3d_slices) { 1324b3e95d5SJeremy L Thompson std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 133dc007f05SJeremy L Thompson std::string P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q"; 1344b3e95d5SJeremy L Thompson std::string option_name = (is_input ? "inputs" : "outputs"); 1354b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 1364b3e95d5SJeremy L Thompson CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 1374b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 1384b3e95d5SJeremy L Thompson CeedBasis_Cuda_shared *basis_data; 1394b3e95d5SJeremy L Thompson CeedBasis basis; 1404b3e95d5SJeremy L Thompson 1414b3e95d5SJeremy L Thompson code << " // -- " << (is_input ? "Input" : "Output") << " field " << i << "\n"; 1424b3e95d5SJeremy L Thompson 1434b3e95d5SJeremy L Thompson // Get field data 1444b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 1454b3e95d5SJeremy L Thompson if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 1464b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1474b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 1484b3e95d5SJeremy L Thompson } 149681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 1504b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 1514b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 1524b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 153dc007f05SJeremy L Thompson if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 154dc007f05SJeremy L Thompson else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d)); 1554b3e95d5SJeremy L Thompson } 1564b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 1574b3e95d5SJeremy L Thompson 1584b3e95d5SJeremy L Thompson // Set field constants 1594b3e95d5SJeremy L Thompson if (eval_mode != CEED_EVAL_WEIGHT) { 1604b3e95d5SJeremy L Thompson code << " const CeedInt " << P_name << " = " << (basis == CEED_BASIS_NONE ? Q_1d : P_1d) << ";\n"; 1614b3e95d5SJeremy L Thompson code << " const CeedInt num_comp" << var_suffix << " = " << num_comp << ";\n"; 1624b3e95d5SJeremy L Thompson } 1634b3e95d5SJeremy L Thompson 1644b3e95d5SJeremy L Thompson // Load basis data 1654b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 1664b3e95d5SJeremy L Thompson switch (eval_mode) { 1674b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 1684b3e95d5SJeremy L Thompson break; 1694b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 1708b97b69aSJeremy L Thompson if (is_at_points) { 1718b97b69aSJeremy L Thompson // AtPoints 1728b97b69aSJeremy L Thompson if (!basis_data->d_chebyshev_interp_1d) { 1738b97b69aSJeremy L Thompson CeedSize interp_bytes; 1748b97b69aSJeremy L Thompson CeedScalar *chebyshev_interp_1d; 1758b97b69aSJeremy L Thompson 1768b97b69aSJeremy L Thompson interp_bytes = P_1d * Q_1d * sizeof(CeedScalar); 1778b97b69aSJeremy L Thompson CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 1788b97b69aSJeremy L Thompson CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 1798b97b69aSJeremy L Thompson CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes)); 1808b97b69aSJeremy L Thompson CeedCallCuda(CeedBasisReturnCeed(basis), 1818b97b69aSJeremy L Thompson cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 1828b97b69aSJeremy L Thompson CeedCallBackend(CeedFree(&chebyshev_interp_1d)); 1838b97b69aSJeremy L Thompson } 1848b97b69aSJeremy L Thompson if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d; 1858b97b69aSJeremy L Thompson else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d; 1868b97b69aSJeremy L Thompson } else { 1878b97b69aSJeremy L Thompson // Standard quadrature 1884b3e95d5SJeremy L Thompson if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 1894b3e95d5SJeremy L Thompson else data->B.outputs[i] = basis_data->d_interp_1d; 1908b97b69aSJeremy L Thompson } 191dc007f05SJeremy L Thompson code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n"; 192f815fac9SJeremy L Thompson code << " LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 1934b3e95d5SJeremy L Thompson break; 1944b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 1958b97b69aSJeremy L Thompson if (is_at_points) { 1968b97b69aSJeremy L Thompson // AtPoints 1978b97b69aSJeremy L Thompson if (!basis_data->d_chebyshev_interp_1d) { 1988b97b69aSJeremy L Thompson CeedSize interp_bytes; 1998b97b69aSJeremy L Thompson CeedScalar *chebyshev_interp_1d; 2008b97b69aSJeremy L Thompson 2018b97b69aSJeremy L Thompson interp_bytes = P_1d * Q_1d * sizeof(CeedScalar); 2028b97b69aSJeremy L Thompson CeedCallBackend(CeedCalloc(P_1d * Q_1d, &chebyshev_interp_1d)); 2038b97b69aSJeremy L Thompson CeedCallBackend(CeedBasisGetChebyshevInterp1D(basis, chebyshev_interp_1d)); 2048b97b69aSJeremy L Thompson CeedCallCuda(CeedBasisReturnCeed(basis), cudaMalloc((void **)&basis_data->d_chebyshev_interp_1d, interp_bytes)); 2058b97b69aSJeremy L Thompson CeedCallCuda(CeedBasisReturnCeed(basis), 2068b97b69aSJeremy L Thompson cudaMemcpy(basis_data->d_chebyshev_interp_1d, chebyshev_interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 2078b97b69aSJeremy L Thompson CeedCallBackend(CeedFree(&chebyshev_interp_1d)); 2088b97b69aSJeremy L Thompson } 2098b97b69aSJeremy L Thompson if (is_input) data->B.inputs[i] = basis_data->d_chebyshev_interp_1d; 2108b97b69aSJeremy L Thompson else data->B.outputs[i] = basis_data->d_chebyshev_interp_1d; 2118b97b69aSJeremy L Thompson } else { 2128b97b69aSJeremy L Thompson // Standard quadrature 2134b3e95d5SJeremy L Thompson if (is_input) data->B.inputs[i] = basis_data->d_interp_1d; 2144b3e95d5SJeremy L Thompson else data->B.outputs[i] = basis_data->d_interp_1d; 2158b97b69aSJeremy L Thompson } 216dc007f05SJeremy L Thompson if (is_tensor) { 217dc007f05SJeremy L Thompson code << " __shared__ CeedScalar s_B" << var_suffix << "[" << P_name << "*" << Q_name << "];\n"; 218f815fac9SJeremy L Thompson code << " LoadMatrix<" << P_name << ", " << Q_name << ">(data, B." << option_name << "[" << i << "], s_B" << var_suffix << ");\n"; 219dc007f05SJeremy L Thompson } 2208b97b69aSJeremy L Thompson if (is_at_points) break; // No G mat for AtPoints 2214b3e95d5SJeremy L Thompson if (use_3d_slices) { 2224b3e95d5SJeremy L Thompson if (is_input) data->G.inputs[i] = basis_data->d_collo_grad_1d; 2234b3e95d5SJeremy L Thompson else data->G.outputs[i] = basis_data->d_collo_grad_1d; 224dc007f05SJeremy L Thompson code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n"; 225f815fac9SJeremy L Thompson code << " LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 2264b3e95d5SJeremy L Thompson } else { 2274b3e95d5SJeremy L Thompson bool has_collo_grad = basis_data->d_collo_grad_1d; 2284b3e95d5SJeremy L Thompson 2294b3e95d5SJeremy L Thompson if (is_input) data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 2304b3e95d5SJeremy L Thompson else data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 2314b3e95d5SJeremy L Thompson if (has_collo_grad) { 232dc007f05SJeremy L Thompson code << " __shared__ CeedScalar s_G" << var_suffix << "[" << Q_name << "*" << Q_name << "];\n"; 233f815fac9SJeremy L Thompson code << " LoadMatrix<" << Q_name << ", " << Q_name << ">(data, G." << option_name << "[" << i << "], s_G" << var_suffix << ");\n"; 2344b3e95d5SJeremy L Thompson } else { 235dc007f05SJeremy L Thompson code << " __shared__ CeedScalar s_G" << var_suffix << "[" << P_name << "*" << Q_name << (is_tensor ? "" : "*dim") << "];\n"; 236dc007f05SJeremy L Thompson code << " LoadMatrix<" << P_name << ", " << Q_name << (is_tensor ? "" : "*dim") << ">(data, G." << option_name << "[" << i << "], s_G" 237dc007f05SJeremy L Thompson << var_suffix << ");\n"; 2384b3e95d5SJeremy L Thompson } 2394b3e95d5SJeremy L Thompson } 2404b3e95d5SJeremy L Thompson break; 2414b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 2424b3e95d5SJeremy L Thompson break; // No action 2434b3e95d5SJeremy L Thompson // LCOV_EXCL_START 2444b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 2454b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 2464b3e95d5SJeremy L Thompson break; // TODO: Not implemented 2474b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 2484b3e95d5SJeremy L Thompson } 2498b97b69aSJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 2504b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 2514b3e95d5SJeremy L Thompson } 2524b3e95d5SJeremy L Thompson 2534b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 2544b3e95d5SJeremy L Thompson // Restriction 2554b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 2564b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelRestriction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim, 257e93651e5SJeremy L Thompson CeedInt field_input_buffer[], CeedOperatorField op_field, CeedQFunctionField qf_field, 258dc007f05SJeremy L Thompson CeedInt Q_1d, bool is_input, bool is_tensor, bool is_at_points, bool use_3d_slices) { 2594b3e95d5SJeremy L Thompson std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 260dc007f05SJeremy L Thompson std::string P_name = (is_tensor ? "P_1d" : "P") + var_suffix; 2614b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 2624b3e95d5SJeremy L Thompson CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 2634b3e95d5SJeremy L Thompson CeedSize l_size; 264f815fac9SJeremy L Thompson CeedRestrictionType rstr_type = CEED_RESTRICTION_STANDARD; 2654b3e95d5SJeremy L Thompson CeedElemRestriction_Cuda *rstr_data; 2664b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 2674b3e95d5SJeremy L Thompson CeedBasis basis; 2684b3e95d5SJeremy L Thompson 2694b3e95d5SJeremy L Thompson // Get field data 2704b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 2714b3e95d5SJeremy L Thompson if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 272f815fac9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type)); 2734b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 2744b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 2754b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 2764b3e95d5SJeremy L Thompson } 2774b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 2784b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 279dc007f05SJeremy L Thompson if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 280dc007f05SJeremy L Thompson else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d)); 2814b3e95d5SJeremy L Thompson } 282681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 2834b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 2844b3e95d5SJeremy L Thompson 2854b3e95d5SJeremy L Thompson // Restriction 2864b3e95d5SJeremy L Thompson if (is_input) { 2874b3e95d5SJeremy L Thompson // Input 288e93651e5SJeremy L Thompson if (field_input_buffer[i] != i) { 289e93651e5SJeremy L Thompson std::string buffer_name = "r_e_in_" + std::to_string(field_input_buffer[i]); 290e93651e5SJeremy L Thompson 291e93651e5SJeremy L Thompson // Restriction was already done for previous input 292e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = " << buffer_name << ";\n"; 2938b97b69aSJeremy L Thompson } else if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_3d_slices && is_at_points)) { 2948b97b69aSJeremy L Thompson if (eval_mode == CEED_EVAL_NONE && rstr_type != CEED_RESTRICTION_POINTS) { 295e93651e5SJeremy L Thompson // No basis action, so r_e_in_* in also r_q_in_* and needs to be allocated 2964b3e95d5SJeremy L Thompson code << " CeedScalar r_e" << var_suffix << "[num_comp" << var_suffix << "*" << P_name << "];\n"; 2978b97b69aSJeremy L Thompson } else if (rstr_type != CEED_RESTRICTION_POINTS) { 298e93651e5SJeremy L Thompson // Otherwise we're using the scratch space 299e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 300e93651e5SJeremy L Thompson } 301f815fac9SJeremy L Thompson switch (rstr_type) { 302f815fac9SJeremy L Thompson case CEED_RESTRICTION_STANDARD: { 3034b3e95d5SJeremy L Thompson CeedInt comp_stride; 3044b3e95d5SJeremy L Thompson 3054b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 3064b3e95d5SJeremy L Thompson code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 3074b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 3084b3e95d5SJeremy L Thompson code << " // CompStride: " << comp_stride << "\n"; 3094b3e95d5SJeremy L Thompson data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 310dc007f05SJeremy L Thompson code << " ReadLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name 311dc007f05SJeremy L Thompson << ">(data, l_size" << var_suffix << ", elem, indices.inputs[" << i << "], d" << var_suffix << ", r_e" << var_suffix << ");\n"; 312f815fac9SJeremy L Thompson break; 313f815fac9SJeremy L Thompson } 314f815fac9SJeremy L Thompson case CEED_RESTRICTION_STRIDED: { 3154b3e95d5SJeremy L Thompson bool has_backend_strides; 3164b3e95d5SJeremy L Thompson CeedInt num_elem; 3174b3e95d5SJeremy L Thompson 3184b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 3194b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 3204b3e95d5SJeremy L Thompson CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 3214b3e95d5SJeremy L Thompson 3224b3e95d5SJeremy L Thompson if (!has_backend_strides) { 3234b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 3244b3e95d5SJeremy L Thompson } 3254b3e95d5SJeremy L Thompson code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 326dc007f05SJeremy L Thompson code << " ReadLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", " 327dc007f05SJeremy L Thompson << strides[1] << ", " << strides[2] << ">(data, elem, d" << var_suffix << ", r_e" << var_suffix << ");\n"; 328f815fac9SJeremy L Thompson break; 329f815fac9SJeremy L Thompson } 3308b97b69aSJeremy L Thompson case CEED_RESTRICTION_POINTS: { 3318b97b69aSJeremy L Thompson CeedInt comp_stride; 3328b97b69aSJeremy L Thompson 3338b97b69aSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 3348b97b69aSJeremy L Thompson code << " const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n"; 3358b97b69aSJeremy L Thompson data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 3368b97b69aSJeremy L Thompson break; 3378b97b69aSJeremy L Thompson } 338f815fac9SJeremy L Thompson // LCOV_EXCL_START 339f815fac9SJeremy L Thompson case CEED_RESTRICTION_ORIENTED: 340f815fac9SJeremy L Thompson case CEED_RESTRICTION_CURL_ORIENTED: 341f815fac9SJeremy L Thompson break; // TODO: Not implemented 342f815fac9SJeremy L Thompson // LCOV_EXCL_STOP 3434b3e95d5SJeremy L Thompson } 3444b3e95d5SJeremy L Thompson } 3454b3e95d5SJeremy L Thompson } else { 3464b3e95d5SJeremy L Thompson // Output 347f815fac9SJeremy L Thompson switch (rstr_type) { 348f815fac9SJeremy L Thompson case CEED_RESTRICTION_STANDARD: { 3494b3e95d5SJeremy L Thompson CeedInt comp_stride; 3504b3e95d5SJeremy L Thompson 3514b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 3524b3e95d5SJeremy L Thompson code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 3534b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 3544b3e95d5SJeremy L Thompson code << " // CompStride: " << comp_stride << "\n"; 3554b3e95d5SJeremy L Thompson data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets; 356dc007f05SJeremy L Thompson code << " WriteLVecStandard" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << comp_stride << ", " << P_name 357dc007f05SJeremy L Thompson << ">(data, l_size" << var_suffix << ", elem, indices.outputs[" << i << "], r_e" << var_suffix << ", d" << var_suffix << ");\n"; 358f815fac9SJeremy L Thompson break; 359f815fac9SJeremy L Thompson } 360f815fac9SJeremy L Thompson case CEED_RESTRICTION_STRIDED: { 3614b3e95d5SJeremy L Thompson bool has_backend_strides; 3624b3e95d5SJeremy L Thompson CeedInt num_elem; 3634b3e95d5SJeremy L Thompson 3644b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 3654b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 3664b3e95d5SJeremy L Thompson CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 3674b3e95d5SJeremy L Thompson 3684b3e95d5SJeremy L Thompson if (!has_backend_strides) { 3694b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 3704b3e95d5SJeremy L Thompson } 3714b3e95d5SJeremy L Thompson code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 372dc007f05SJeremy L Thompson code << " WriteLVecStrided" << (is_tensor ? dim : 1) << "d<num_comp" << var_suffix << ", " << P_name << ", " << strides[0] << ", " 373dc007f05SJeremy L Thompson << strides[1] << ", " << strides[2] << ">(data, elem, r_e" << var_suffix << ", d" << var_suffix << ");\n"; 374f815fac9SJeremy L Thompson break; 375f815fac9SJeremy L Thompson } 3768b97b69aSJeremy L Thompson case CEED_RESTRICTION_POINTS: 3778b97b69aSJeremy L Thompson data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets; 3788b97b69aSJeremy L Thompson break; 379f815fac9SJeremy L Thompson // LCOV_EXCL_START 380f815fac9SJeremy L Thompson case CEED_RESTRICTION_ORIENTED: 381f815fac9SJeremy L Thompson case CEED_RESTRICTION_CURL_ORIENTED: 382f815fac9SJeremy L Thompson break; // TODO: Not implemented 383f815fac9SJeremy L Thompson // LCOV_EXCL_STOP 3844b3e95d5SJeremy L Thompson } 3854b3e95d5SJeremy L Thompson } 386681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 3874b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 3884b3e95d5SJeremy L Thompson } 3894b3e95d5SJeremy L Thompson 3904b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 3914b3e95d5SJeremy L Thompson // Basis 3924b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 3934b3e95d5SJeremy L Thompson static int CeedOperatorBuildKernelBasis_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt i, CeedInt dim, 394dc007f05SJeremy L Thompson CeedOperatorField op_field, CeedQFunctionField qf_field, CeedInt Q_1d, bool is_input, bool is_tensor, 3958b97b69aSJeremy L Thompson bool is_at_points, bool use_3d_slices) { 3964b3e95d5SJeremy L Thompson std::string var_suffix = (is_input ? "_in_" : "_out_") + std::to_string(i); 397dc007f05SJeremy L Thompson std::string P_name = (is_tensor ? "P_1d" : "P") + var_suffix, Q_name = is_tensor ? "Q_1d" : "Q"; 3984b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 3994b3e95d5SJeremy L Thompson CeedInt elem_size = 0, num_comp = 0, P_1d = 0; 4004b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 4014b3e95d5SJeremy L Thompson CeedBasis basis; 4024b3e95d5SJeremy L Thompson 4034b3e95d5SJeremy L Thompson // Get field data 4044b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_field, &elem_rstr)); 4054b3e95d5SJeremy L Thompson if (elem_rstr != CEED_ELEMRESTRICTION_NONE) { 4064b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 4074b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 4084b3e95d5SJeremy L Thompson } 409681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 4104b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetBasis(op_field, &basis)); 4114b3e95d5SJeremy L Thompson if (basis != CEED_BASIS_NONE) { 412dc007f05SJeremy L Thompson if (is_tensor) CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 413dc007f05SJeremy L Thompson else CeedCallBackend(CeedBasisGetNumNodes(basis, &P_1d)); 4144b3e95d5SJeremy L Thompson } 4154b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_field, &eval_mode)); 4164b3e95d5SJeremy L Thompson 4174b3e95d5SJeremy L Thompson // Basis 4184b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 4194b3e95d5SJeremy L Thompson if (is_input) { 4204b3e95d5SJeremy L Thompson switch (eval_mode) { 4214b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 4228b97b69aSJeremy L Thompson if (!use_3d_slices && !is_at_points) { 4234b3e95d5SJeremy L Thompson code << " CeedScalar *r_q" << var_suffix << " = r_e" << var_suffix << ";\n"; 4244b3e95d5SJeremy L Thompson } 4254b3e95d5SJeremy L Thompson break; 4264b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 4278b97b69aSJeremy L Thompson if (is_at_points) { 428dc007f05SJeremy L Thompson std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d"; 429dc007f05SJeremy L Thompson 430dc007f05SJeremy L Thompson code << " CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n"; 431dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" 432dc007f05SJeremy L Thompson << var_suffix << ", r_c" << var_suffix << ");\n"; 4338b97b69aSJeremy L Thompson } else { 434dc007f05SJeremy L Thompson std::string function_name = is_tensor ? ((dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d") : "InterpNonTensor"; 435dc007f05SJeremy L Thompson 436dc007f05SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n"; 437dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" 438dc007f05SJeremy L Thompson << var_suffix << ", r_q" << var_suffix << ");\n"; 4398b97b69aSJeremy L Thompson } 4404b3e95d5SJeremy L Thompson break; 4414b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 4428b97b69aSJeremy L Thompson if (is_at_points) { 443dc007f05SJeremy L Thompson std::string function_name = (dim == 1 ? "Interp" : "InterpTensor") + std::to_string(dim) + "d"; 444dc007f05SJeremy L Thompson 445dc007f05SJeremy L Thompson code << " CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n"; 446dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" 447dc007f05SJeremy L Thompson << var_suffix << ", r_c" << var_suffix << ");\n"; 4488b97b69aSJeremy L Thompson } else if (use_3d_slices) { 449dc007f05SJeremy L Thompson std::string function_name = (dim > 1 ? "InterpTensor" : "Interp") + std::to_string(dim) + "d"; 450dc007f05SJeremy L Thompson 4514b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 452dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" 453dc007f05SJeremy L Thompson << var_suffix << ", r_q" << var_suffix << ");\n"; 454dc007f05SJeremy L Thompson } else if (is_tensor) { 455dc007f05SJeremy L Thompson bool is_collocated = dim == 3 && Q_1d >= P_1d; 456dc007f05SJeremy L Thompson std::string function_name = (dim == 1 ? "Grad" : (is_collocated ? "GradTensorCollocated" : "GradTensor")) + std::to_string(dim) + "d"; 457dc007f05SJeremy L Thompson 458dc007f05SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (dim >= 3 ? Q_name : "1") << "];\n"; 459dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix << ", s_B" 460dc007f05SJeremy L Thompson << var_suffix << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n"; 4614b3e95d5SJeremy L Thompson } else { 462dc007f05SJeremy L Thompson std::string function_name = "GradNonTensor"; 463dc007f05SJeremy L Thompson 464dc007f05SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 465dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_e" << var_suffix 466dc007f05SJeremy L Thompson << ", s_G" << var_suffix << ", r_q" << var_suffix << ");\n"; 4674b3e95d5SJeremy L Thompson } 4684b3e95d5SJeremy L Thompson break; 4694b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: { 4708b97b69aSJeremy L Thompson if (is_at_points) { 4718b97b69aSJeremy L Thompson code << " // Nothing to do AtPoints\n"; 4728b97b69aSJeremy L Thompson } else { 4734b3e95d5SJeremy L Thompson CeedBasis_Cuda_shared *basis_data; 474dc007f05SJeremy L Thompson std::string function_name = is_tensor ? ((dim == 1 ? "Weight" : "WeightTensor") + std::to_string(dim) + "d") : "WeightNonTensor"; 4754b3e95d5SJeremy L Thompson 476dc007f05SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n"; 4774b3e95d5SJeremy L Thompson CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 4784b3e95d5SJeremy L Thompson data->W = basis_data->d_q_weight_1d; 479dc007f05SJeremy L Thompson code << " " << function_name << "<" << Q_name << ">(data, W, r_q" << var_suffix << ");\n"; 4808b97b69aSJeremy L Thompson } 4814b3e95d5SJeremy L Thompson break; 4824b3e95d5SJeremy L Thompson } 4834b3e95d5SJeremy L Thompson // LCOV_EXCL_START 4844b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 4854b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 4864b3e95d5SJeremy L Thompson break; // TODO: Not implemented 4874b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 4884b3e95d5SJeremy L Thompson } 4894b3e95d5SJeremy L Thompson } else { 4904b3e95d5SJeremy L Thompson switch (eval_mode) { 4914b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 4924b3e95d5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_q" << var_suffix << ";\n"; 4934b3e95d5SJeremy L Thompson break; // No action 4944b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 495e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 4968b97b69aSJeremy L Thompson if (is_at_points) { 497dc007f05SJeremy L Thompson std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d"; 498dc007f05SJeremy L Thompson 499dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B" 500dc007f05SJeremy L Thompson << var_suffix << ", r_e" << var_suffix << ");\n"; 5018b97b69aSJeremy L Thompson } else { 502dc007f05SJeremy L Thompson std::string function_name = 503dc007f05SJeremy L Thompson is_tensor ? ((dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d") : "InterpTransposeNonTensor"; 504dc007f05SJeremy L Thompson 505dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B" 506dc007f05SJeremy L Thompson << var_suffix << ", r_e" << var_suffix << ");\n"; 5078b97b69aSJeremy L Thompson } 5084b3e95d5SJeremy L Thompson break; 5094b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 510e93651e5SJeremy L Thompson code << " CeedScalar *r_e" << var_suffix << " = r_e_scratch;\n"; 5118b97b69aSJeremy L Thompson if (is_at_points) { 512dc007f05SJeremy L Thompson std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d"; 513dc007f05SJeremy L Thompson 514dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_c" << var_suffix << ", s_B" 515dc007f05SJeremy L Thompson << var_suffix << ", r_e" << var_suffix << ");\n"; 5168b97b69aSJeremy L Thompson } else if (use_3d_slices) { 517dc007f05SJeremy L Thompson std::string function_name = (dim == 1 ? "InterpTranspose" : "InterpTransposeTensor") + std::to_string(dim) + "d"; 518dc007f05SJeremy L Thompson 519dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B" 520dc007f05SJeremy L Thompson << var_suffix << ", r_e" << var_suffix << ");\n"; 521dc007f05SJeremy L Thompson } else if (is_tensor) { 522dc007f05SJeremy L Thompson bool is_collocated = dim == 3 && Q_1d >= P_1d; 523dc007f05SJeremy L Thompson std::string function_name = 524dc007f05SJeremy L Thompson (dim == 1 ? "GradTranspose" : (is_collocated ? "GradTransposeTensorCollocated" : "GradTransposeTensor")) + std::to_string(dim) + "d"; 525dc007f05SJeremy L Thompson 526dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix << ", s_B" 527dc007f05SJeremy L Thompson << var_suffix << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n"; 5284b3e95d5SJeremy L Thompson } else { 529dc007f05SJeremy L Thompson std::string function_name = "GradTransposeNonTensor"; 530dc007f05SJeremy L Thompson 531dc007f05SJeremy L Thompson code << " " << function_name << "<num_comp" << var_suffix << ", dim, " << P_name << ", " << Q_name << ">(data, r_q" << var_suffix 532dc007f05SJeremy L Thompson << ", s_G" << var_suffix << ", r_e" << var_suffix << ");\n"; 5334b3e95d5SJeremy L Thompson } 5344b3e95d5SJeremy L Thompson break; 5354b3e95d5SJeremy L Thompson // LCOV_EXCL_START 5364b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 5374b3e95d5SJeremy L Thompson break; // Should not occur 5384b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 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 } 544681d0ea7SJeremy L Thompson CeedCallBackend(CeedBasisDestroy(&basis)); 5454b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 5464b3e95d5SJeremy L Thompson } 5474b3e95d5SJeremy L Thompson 5484b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 5494b3e95d5SJeremy L Thompson // QFunction 5504b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 5518b97b69aSJeremy L Thompson static int CeedOperatorBuildKernelQFunction_Cuda_gen(std::ostringstream &code, CeedOperator_Cuda_gen *data, CeedInt dim, CeedInt max_num_points, 5528b97b69aSJeremy L Thompson CeedInt num_input_fields, CeedOperatorField *op_input_fields, 5538b97b69aSJeremy L Thompson CeedQFunctionField *qf_input_fields, CeedInt num_output_fields, 5548b97b69aSJeremy L Thompson CeedOperatorField *op_output_fields, CeedQFunctionField *qf_output_fields, 555dc007f05SJeremy L Thompson std::string qfunction_name, CeedInt Q_1d, bool is_tensor, bool is_at_points, 556dc007f05SJeremy L Thompson bool use_3d_slices) { 557dc007f05SJeremy L Thompson std::string Q_name = is_tensor ? "Q_1d" : "Q"; 5584b3e95d5SJeremy L Thompson CeedEvalMode eval_mode = CEED_EVAL_NONE; 5594b3e95d5SJeremy L Thompson CeedElemRestriction elem_rstr; 5604b3e95d5SJeremy L Thompson 5618b97b69aSJeremy L Thompson // Setup output arrays 5624b3e95d5SJeremy L Thompson code << "\n // -- Output field setup\n"; 5634b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 5644b3e95d5SJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 5654b3e95d5SJeremy L Thompson 5664b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 5674b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 5688b97b69aSJeremy L Thompson switch (eval_mode) { 5698b97b69aSJeremy L Thompson case CEED_EVAL_NONE: 5708b97b69aSJeremy L Thompson if (is_at_points) { 5718b97b69aSJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "];\n"; 5728b97b69aSJeremy L Thompson } else { 573dc007f05SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n"; 5744b3e95d5SJeremy L Thompson } 5758b97b69aSJeremy L Thompson break; 5768b97b69aSJeremy L Thompson case CEED_EVAL_INTERP: 5778b97b69aSJeremy L Thompson if (is_at_points) { 5788b97b69aSJeremy L Thompson // Accumulator for point data 579dc007f05SJeremy L Thompson code << " CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "];\n"; 580dc007f05SJeremy L Thompson code << " for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n"; 5818b97b69aSJeremy L Thompson code << " r_c" << var_suffix << "[i] = 0.0;\n"; 5828b97b69aSJeremy L Thompson code << " }\n"; 5838b97b69aSJeremy L Thompson } else { 584dc007f05SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n"; 5858b97b69aSJeremy L Thompson } 5868b97b69aSJeremy L Thompson break; 5878b97b69aSJeremy L Thompson case CEED_EVAL_GRAD: 5888b97b69aSJeremy L Thompson if (is_at_points) { 5898b97b69aSJeremy L Thompson // Accumulator for point data 590dc007f05SJeremy L Thompson code << " CeedScalar r_c" << var_suffix << "[num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "*dim];\n"; 591dc007f05SJeremy L Thompson code << " for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << (dim >= 3 ? Q_name : "1") << "; i++) {\n"; 5928b97b69aSJeremy L Thompson code << " r_c" << var_suffix << "[i] = 0.0;\n"; 5938b97b69aSJeremy L Thompson code << " }\n"; 5948b97b69aSJeremy L Thompson } else if (use_3d_slices) { 5954b3e95d5SJeremy L Thompson // Accumulator for gradient slices 5964b3e95d5SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*" << Q_name << "];\n"; 5974b3e95d5SJeremy L Thompson code << " for (CeedInt i = 0; i < num_comp" << var_suffix << "*" << Q_name << "; i++) {\n"; 5984b3e95d5SJeremy L Thompson code << " r_q" << var_suffix << "[i] = 0.0;\n"; 5994b3e95d5SJeremy L Thompson code << " }\n"; 6004b3e95d5SJeremy L Thompson } else { 601dc007f05SJeremy L Thompson code << " CeedScalar r_q" << var_suffix << "[num_comp" << var_suffix << "*dim*" << (is_tensor && (dim >= 3) ? Q_name : "1") << "];\n"; 6024b3e95d5SJeremy L Thompson } 6038b97b69aSJeremy L Thompson break; 6048b97b69aSJeremy L Thompson case CEED_EVAL_WEIGHT: 6058b97b69aSJeremy L Thompson break; 6068b97b69aSJeremy L Thompson // LCOV_EXCL_START 6078b97b69aSJeremy L Thompson case CEED_EVAL_DIV: 6088b97b69aSJeremy L Thompson case CEED_EVAL_CURL: 6098b97b69aSJeremy L Thompson break; // TODO: Not implemented 6108b97b69aSJeremy L Thompson // LCOV_EXCL_STOP 6114b3e95d5SJeremy L Thompson } 6124b3e95d5SJeremy L Thompson } 6134b3e95d5SJeremy L Thompson 6148b97b69aSJeremy L Thompson if (is_at_points) { 6158b97b69aSJeremy L Thompson // We need to handle batches of points 6168b97b69aSJeremy L Thompson code << "\n // Note: Using batches of points\n"; 6178b97b69aSJeremy L Thompson code << " const CeedInt point_loop_bound = (blockDim.x * blockDim.y) * ceil(1.0 * max_num_points / (blockDim.x * blockDim.y));\n\n"; 6188b97b69aSJeremy L Thompson code << " #pragma unroll\n"; 6198b97b69aSJeremy L Thompson code << " for (CeedInt i = threadIdx.x + threadIdx.y * blockDim.x; i < point_loop_bound; i += blockDim.x * blockDim.y) {\n"; 6208b97b69aSJeremy L Thompson code << " const CeedInt p = i % max_num_points;\n\n"; 6218b97b69aSJeremy L Thompson 6228b97b69aSJeremy L Thompson code << " // -- Coordinates\n"; 6238b97b69aSJeremy L Thompson code << " CeedScalar r_x[dim];\n"; 6248b97b69aSJeremy L Thompson code << " ReadPoint<dim, coords_comp_stride, max_num_points>(data, elem, p, max_num_points, points.indices, points.coords, r_x);\n\n"; 6258b97b69aSJeremy L Thompson 6268b97b69aSJeremy L Thompson code << " // -- Input fields\n"; 6278b97b69aSJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 6288b97b69aSJeremy L Thompson std::string var_suffix = "_in_" + std::to_string(i); 6298b97b69aSJeremy L Thompson 6308b97b69aSJeremy L Thompson code << " // ---- Input field " << i << "\n"; 6318b97b69aSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 6328b97b69aSJeremy L Thompson // Basis action 6338b97b69aSJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 6348b97b69aSJeremy L Thompson switch (eval_mode) { 6358b97b69aSJeremy L Thompson case CEED_EVAL_NONE: 6368b97b69aSJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 6378b97b69aSJeremy L Thompson code << " ReadPoint<num_comp" << var_suffix << ", comp_stride" << var_suffix 6388b97b69aSJeremy L Thompson << ", max_num_points>(data, elem, p, max_num_points, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n"; 6398b97b69aSJeremy L Thompson break; 6408b97b69aSJeremy L Thompson case CEED_EVAL_INTERP: 6418b97b69aSJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 6428b97b69aSJeremy L Thompson code << " InterpAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix 6438b97b69aSJeremy L Thompson << ", r_x, r_s" << var_suffix << ");\n"; 6448b97b69aSJeremy L Thompson break; 6458b97b69aSJeremy L Thompson case CEED_EVAL_GRAD: 6468b97b69aSJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 6478b97b69aSJeremy L Thompson code << " GradAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_c" << var_suffix 6488b97b69aSJeremy L Thompson << ", r_x, r_s" << var_suffix << ");\n"; 6498b97b69aSJeremy L Thompson break; 6508b97b69aSJeremy L Thompson case CEED_EVAL_WEIGHT: 6518b97b69aSJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[1];\n"; 6528b97b69aSJeremy L Thompson code << " r_s" << var_suffix << "[0] = 1.0;\n"; 6538b97b69aSJeremy L Thompson break; 6548b97b69aSJeremy L Thompson // LCOV_EXCL_START 6558b97b69aSJeremy L Thompson case CEED_EVAL_DIV: 6568b97b69aSJeremy L Thompson case CEED_EVAL_CURL: 6578b97b69aSJeremy L Thompson break; // TODO: Not implemented 6588b97b69aSJeremy L Thompson // LCOV_EXCL_STOP 6598b97b69aSJeremy L Thompson } 6608b97b69aSJeremy L Thompson } 6618b97b69aSJeremy L Thompson code << "\n // -- Output fields\n"; 6628b97b69aSJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 6638b97b69aSJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 6648b97b69aSJeremy L Thompson 6658b97b69aSJeremy L Thompson code << " // ---- Output field " << i << "\n"; 6668b97b69aSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 6678b97b69aSJeremy L Thompson // Basis action 6688b97b69aSJeremy L Thompson switch (eval_mode) { 6698b97b69aSJeremy L Thompson case CEED_EVAL_NONE: 6708b97b69aSJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 6718b97b69aSJeremy L Thompson break; 6728b97b69aSJeremy L Thompson case CEED_EVAL_INTERP: 6738b97b69aSJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 6748b97b69aSJeremy L Thompson break; 6758b97b69aSJeremy L Thompson case CEED_EVAL_GRAD: 6768b97b69aSJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 6778b97b69aSJeremy L Thompson break; 6788b97b69aSJeremy L Thompson // LCOV_EXCL_START 6798b97b69aSJeremy L Thompson case CEED_EVAL_WEIGHT: 6808b97b69aSJeremy L Thompson break; // Should not occur 6818b97b69aSJeremy L Thompson case CEED_EVAL_DIV: 6828b97b69aSJeremy L Thompson case CEED_EVAL_CURL: 6838b97b69aSJeremy L Thompson break; // TODO: Not implemented 6848b97b69aSJeremy L Thompson // LCOV_EXCL_STOP 6858b97b69aSJeremy L Thompson } 6868b97b69aSJeremy L Thompson } 6878b97b69aSJeremy L Thompson 6888b97b69aSJeremy L Thompson } else if (use_3d_slices) { 6894b3e95d5SJeremy L Thompson // We treat quadrature points per slice in 3d to save registers 6904b3e95d5SJeremy L Thompson code << "\n // Note: Using planes of 3D elements\n"; 6914b3e95d5SJeremy L Thompson code << " #pragma unroll\n"; 6924b3e95d5SJeremy L Thompson code << " for (CeedInt q = 0; q < " << Q_name << "; q++) {\n"; 6934b3e95d5SJeremy L Thompson code << " // -- Input fields\n"; 6944b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 6954b3e95d5SJeremy L Thompson std::string var_suffix = "_in_" + std::to_string(i); 6964b3e95d5SJeremy L Thompson 6974b3e95d5SJeremy L Thompson code << " // ---- Input field " << i << "\n"; 6984b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 6994b3e95d5SJeremy L Thompson // Basis action 7004b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 7014b3e95d5SJeremy L Thompson switch (eval_mode) { 7024b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 7034b3e95d5SJeremy L Thompson bool is_strided; 7044b3e95d5SJeremy L Thompson 7054b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 7064b3e95d5SJeremy L Thompson 7074b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 7084b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 7094b3e95d5SJeremy L Thompson if (is_strided) { 7104b3e95d5SJeremy L Thompson bool has_backend_strides; 7114b3e95d5SJeremy L Thompson CeedInt num_elem, elem_size; 7124b3e95d5SJeremy L Thompson 7134b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 7144b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 7154b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 7164b3e95d5SJeremy L Thompson CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 7174b3e95d5SJeremy L Thompson 7184b3e95d5SJeremy L Thompson if (!has_backend_strides) { 7194b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 7204b3e95d5SJeremy L Thompson } 7214b3e95d5SJeremy L Thompson code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 722f815fac9SJeremy L Thompson code << " ReadEVecSliceStrided3d<num_comp" << var_suffix << ", " << Q_name << ", " << strides[0] << ", " << strides[1] << ", " 7234b3e95d5SJeremy L Thompson << strides[2] << ">(data, elem, q, d" << var_suffix << ", r_s" << var_suffix << ");\n"; 7244b3e95d5SJeremy L Thompson } else { 7254b3e95d5SJeremy L Thompson CeedSize l_size = 0; 7264b3e95d5SJeremy L Thompson CeedInt comp_stride; 7274b3e95d5SJeremy L Thompson CeedElemRestriction_Cuda *rstr_data; 7284b3e95d5SJeremy L Thompson 7294b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 7304b3e95d5SJeremy L Thompson code << " const CeedInt l_size" << var_suffix << " = " << l_size << ";\n"; 7314b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 7324b3e95d5SJeremy L Thompson code << " // CompStride: " << comp_stride << "\n"; 7334b3e95d5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 7344b3e95d5SJeremy L Thompson data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 735f815fac9SJeremy L Thompson code << " ReadEVecSliceStandard3d<num_comp" << var_suffix << ", " << comp_stride << ", " << Q_name << ">(data, l_size" << var_suffix 7364b3e95d5SJeremy L Thompson << ", elem, q, indices.inputs[" << i << "], d" << var_suffix << ", r_s" << var_suffix << ");\n"; 7374b3e95d5SJeremy L Thompson } 738681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 7394b3e95d5SJeremy L Thompson break; 7404b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 7414b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 7424b3e95d5SJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) {\n"; 7434b3e95d5SJeremy L Thompson code << " r_s" << var_suffix << "[j] = r_q" << var_suffix << "[q + j*" << Q_name << "];\n"; 7444b3e95d5SJeremy L Thompson code << " }\n"; 7454b3e95d5SJeremy L Thompson break; 7464b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 7474b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 748f815fac9SJeremy L Thompson code << " GradColloSlice3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_q" << var_suffix << ", s_G" << var_suffix 749f815fac9SJeremy L Thompson << ", r_s" << var_suffix << ");\n"; 7504b3e95d5SJeremy L Thompson break; 7514b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 7524b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[1];\n"; 7534b3e95d5SJeremy L Thompson code << " r_s" << var_suffix << "[0] = r_q" << var_suffix << "[q];\n"; 7548b97b69aSJeremy L Thompson break; 7554b3e95d5SJeremy L Thompson // LCOV_EXCL_START 7564b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 7574b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 7584b3e95d5SJeremy L Thompson break; // TODO: Not implemented 7594b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 7604b3e95d5SJeremy L Thompson } 7614b3e95d5SJeremy L Thompson } 7624b3e95d5SJeremy L Thompson code << "\n // -- Output fields\n"; 7634b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 7644b3e95d5SJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 7654b3e95d5SJeremy L Thompson 7664b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 7674b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 7684b3e95d5SJeremy L Thompson // Basis action 7694b3e95d5SJeremy L Thompson switch (eval_mode) { 7704b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 7714b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 7728b97b69aSJeremy L Thompson break; 7734b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 7744b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "];\n"; 7754b3e95d5SJeremy L Thompson break; 7764b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 7774b3e95d5SJeremy L Thompson code << " CeedScalar r_s" << var_suffix << "[num_comp" << var_suffix << "*dim];\n"; 7784b3e95d5SJeremy L Thompson break; 7794b3e95d5SJeremy L Thompson // LCOV_EXCL_START 7804b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 7814b3e95d5SJeremy L Thompson break; // Should not occur 7824b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 7834b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 7844b3e95d5SJeremy L Thompson break; // TODO: Not implemented 7854b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 7864b3e95d5SJeremy L Thompson } 7874b3e95d5SJeremy L Thompson } 7884b3e95d5SJeremy L Thompson } else { 7894b3e95d5SJeremy L Thompson code << "\n // Note: Using full elements\n"; 7904b3e95d5SJeremy L Thompson code << " {\n"; 7914b3e95d5SJeremy L Thompson code << " // -- Input fields\n"; 7924b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 7934b3e95d5SJeremy L Thompson code << " // ---- Input field " << i << "\n"; 7944b3e95d5SJeremy L Thompson code << " CeedScalar *r_s_in_" << i << " = r_q_in_" << i << ";\n"; 7954b3e95d5SJeremy L Thompson } 7964b3e95d5SJeremy L Thompson code << " // -- Output fields\n"; 7974b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 7984b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 7994b3e95d5SJeremy L Thompson code << " CeedScalar *r_s_out_" << i << " = r_q_out_" << i << ";\n"; 8004b3e95d5SJeremy L Thompson } 8014b3e95d5SJeremy L Thompson } 8024b3e95d5SJeremy L Thompson 8034b3e95d5SJeremy L Thompson // Input and output buffers 8044b3e95d5SJeremy L Thompson code << "\n // -- QFunction inputs and outputs\n"; 8054b3e95d5SJeremy L Thompson code << " // ---- Inputs\n"; 8064b3e95d5SJeremy L Thompson code << " CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n"; 8074b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 8084b3e95d5SJeremy L Thompson code << " // ------ Input field " << i << "\n"; 8094b3e95d5SJeremy L Thompson code << " inputs[" << i << "] = r_s_in_" << i << ";\n"; 8104b3e95d5SJeremy L Thompson } 8114b3e95d5SJeremy L Thompson code << " // ---- Outputs\n"; 8124b3e95d5SJeremy L Thompson code << " CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n"; 8134b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 8144b3e95d5SJeremy L Thompson code << " // ------ Output field " << i << "\n"; 8154b3e95d5SJeremy L Thompson code << " outputs[" << i << "] = r_s_out_" << i << ";\n"; 8164b3e95d5SJeremy L Thompson } 8174b3e95d5SJeremy L Thompson 8184b3e95d5SJeremy L Thompson // Apply QFunction 8194b3e95d5SJeremy L Thompson code << "\n // -- Apply QFunction\n"; 8204b3e95d5SJeremy L Thompson code << " " << qfunction_name << "(ctx, "; 821dc007f05SJeremy L Thompson if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) { 8224b3e95d5SJeremy L Thompson code << "1"; 8234b3e95d5SJeremy L Thompson } else { 824dc007f05SJeremy L Thompson code << Q_name; 8254b3e95d5SJeremy L Thompson } 8264b3e95d5SJeremy L Thompson code << ", inputs, outputs);\n"; 8274b3e95d5SJeremy L Thompson 8288b97b69aSJeremy L Thompson if (is_at_points) { 8298b97b69aSJeremy L Thompson // Map back to coefficients 8308b97b69aSJeremy L Thompson code << "\n // -- Output fields\n"; 8318b97b69aSJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 8328b97b69aSJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 8338b97b69aSJeremy L Thompson std::string P_name = "P_1d" + var_suffix; 8348b97b69aSJeremy L Thompson 8358b97b69aSJeremy L Thompson code << " // ---- Output field " << i << "\n"; 8368b97b69aSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 8378b97b69aSJeremy L Thompson // Basis action 8388b97b69aSJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 8398b97b69aSJeremy L Thompson switch (eval_mode) { 8408b97b69aSJeremy L Thompson case CEED_EVAL_NONE: { 8418b97b69aSJeremy L Thompson CeedInt comp_stride; 8428b97b69aSJeremy L Thompson CeedElemRestriction elem_rstr; 8438b97b69aSJeremy L Thompson 8448b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 8458b97b69aSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 8468b97b69aSJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 8478b97b69aSJeremy L Thompson code << " const CeedInt comp_stride" << var_suffix << " = " << comp_stride << ";\n"; 8488b97b69aSJeremy L Thompson code << " WritePoint<num_comp" << var_suffix << ", comp_stride" << var_suffix 8498b97b69aSJeremy L Thompson << ", max_num_points>(data, elem, i, points.num_per_elem[elem], indices.outputs[" << i << "]" 8508b97b69aSJeremy L Thompson << ", r_s" << var_suffix << ", d" << var_suffix << ");\n"; 8518b97b69aSJeremy L Thompson break; 8528b97b69aSJeremy L Thompson } 8538b97b69aSJeremy L Thompson case CEED_EVAL_INTERP: 8548b97b69aSJeremy L Thompson code << " if (i >= points.num_per_elem[elem]) {\n"; 8558b97b69aSJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << "; j++) r_s" << var_suffix << "[j] = 0.0;\n"; 8568b97b69aSJeremy L Thompson code << " }\n"; 8578b97b69aSJeremy L Thompson code << " InterpTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s" 8588b97b69aSJeremy L Thompson << var_suffix << ", r_x, r_c" << var_suffix << ");\n"; 8598b97b69aSJeremy L Thompson break; 8608b97b69aSJeremy L Thompson case CEED_EVAL_GRAD: 8618b97b69aSJeremy L Thompson code << " if (i >= points.num_per_elem[elem]) {\n"; 8628b97b69aSJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << "*dim; j++) r_s" << var_suffix << "[j] = 0.0;\n"; 8638b97b69aSJeremy L Thompson code << " }\n"; 8648b97b69aSJeremy L Thompson code << " GradTransposeAtPoints" << dim << "d<num_comp" << var_suffix << ", max_num_points, " << Q_name << ">(data, i, r_s" 8658b97b69aSJeremy L Thompson << var_suffix << ", r_x, r_c" << var_suffix << ");\n"; 8668b97b69aSJeremy L Thompson break; 8678b97b69aSJeremy L Thompson // LCOV_EXCL_START 8688b97b69aSJeremy L Thompson case CEED_EVAL_WEIGHT: 8698b97b69aSJeremy L Thompson break; // Should not occur 8708b97b69aSJeremy L Thompson case CEED_EVAL_DIV: 8718b97b69aSJeremy L Thompson case CEED_EVAL_CURL: 8728b97b69aSJeremy L Thompson break; // TODO: Not implemented 8738b97b69aSJeremy L Thompson // LCOV_EXCL_STOP 8748b97b69aSJeremy L Thompson } 8758b97b69aSJeremy L Thompson } 8768b97b69aSJeremy L Thompson } else if (use_3d_slices) { 8774b3e95d5SJeremy L Thompson // Copy or apply transpose grad, if needed 8788b97b69aSJeremy L Thompson code << "\n // -- Output fields\n"; 8794b3e95d5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 8804b3e95d5SJeremy L Thompson std::string var_suffix = "_out_" + std::to_string(i); 8814b3e95d5SJeremy L Thompson std::string P_name = "P_1d" + var_suffix; 8824b3e95d5SJeremy L Thompson 8834b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 8844b3e95d5SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 8854b3e95d5SJeremy L Thompson // Basis action 8864b3e95d5SJeremy L Thompson code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 8874b3e95d5SJeremy L Thompson switch (eval_mode) { 8884b3e95d5SJeremy L Thompson case CEED_EVAL_NONE: 8894b3e95d5SJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 8904b3e95d5SJeremy L Thompson code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 8914b3e95d5SJeremy L Thompson code << " }\n"; 8928b97b69aSJeremy L Thompson break; 8934b3e95d5SJeremy L Thompson case CEED_EVAL_INTERP: 8944b3e95d5SJeremy L Thompson code << " for (CeedInt j = 0; j < num_comp" << var_suffix << " ; j++) {\n"; 8954b3e95d5SJeremy L Thompson code << " r_q" << var_suffix << "[q + j*" << Q_name << "] = r_s" << var_suffix << "[j];\n"; 8964b3e95d5SJeremy L Thompson code << " }\n"; 8974b3e95d5SJeremy L Thompson break; 8984b3e95d5SJeremy L Thompson case CEED_EVAL_GRAD: 899f815fac9SJeremy L Thompson code << " GradColloSliceTranspose3d<num_comp" << var_suffix << ", " << Q_name << ">(data, q, r_s" << var_suffix << ", s_G" 900f815fac9SJeremy L Thompson << var_suffix << ", r_q" << var_suffix << ");\n"; 9014b3e95d5SJeremy L Thompson break; 9024b3e95d5SJeremy L Thompson // LCOV_EXCL_START 9034b3e95d5SJeremy L Thompson case CEED_EVAL_WEIGHT: 9044b3e95d5SJeremy L Thompson break; // Should not occur 9054b3e95d5SJeremy L Thompson case CEED_EVAL_DIV: 9064b3e95d5SJeremy L Thompson case CEED_EVAL_CURL: 9074b3e95d5SJeremy L Thompson break; // TODO: Not implemented 9084b3e95d5SJeremy L Thompson // LCOV_EXCL_STOP 9094b3e95d5SJeremy L Thompson } 9104b3e95d5SJeremy L Thompson } 9114b3e95d5SJeremy L Thompson } 9124b3e95d5SJeremy L Thompson code << " }\n"; 9134b3e95d5SJeremy L Thompson return CEED_ERROR_SUCCESS; 9144b3e95d5SJeremy L Thompson } 9154b3e95d5SJeremy L Thompson 9164b3e95d5SJeremy L Thompson //------------------------------------------------------------------------------ 917b2165e7aSSebastian Grimberg // Build single operator kernel 918ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 919eb7e6cafSJeremy L Thompson extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op) { 9208b97b69aSJeremy L Thompson bool is_tensor = true, is_at_points = false, use_3d_slices = false; 921241a4b83SYohann Ceed ceed; 9228b97b69aSJeremy L Thompson CeedInt Q_1d, num_input_fields, num_output_fields, dim = 1, max_num_points = 0, coords_comp_stride = 0; 923ca735530SJeremy L Thompson CeedQFunctionField *qf_input_fields, *qf_output_fields; 924ca735530SJeremy L Thompson CeedQFunction_Cuda_gen *qf_data; 925ca735530SJeremy L Thompson CeedQFunction qf; 926ca735530SJeremy L Thompson CeedOperatorField *op_input_fields, *op_output_fields; 927ca735530SJeremy L Thompson CeedOperator_Cuda_gen *data; 9284b3e95d5SJeremy L Thompson std::ostringstream code; 9294b3e95d5SJeremy L Thompson 9304b3e95d5SJeremy L Thompson { 9314b3e95d5SJeremy L Thompson bool is_setup_done; 932ca735530SJeremy L Thompson 933ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 934ca735530SJeremy L Thompson if (is_setup_done) return CEED_ERROR_SUCCESS; 9354b3e95d5SJeremy L Thompson } 936ca735530SJeremy L Thompson 937ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 938ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorGetData(op, &data)); 939ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 940ca735530SJeremy L Thompson CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 941ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 942ca735530SJeremy L Thompson CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 943241a4b83SYohann 9444b3e95d5SJeremy L Thompson // Get operator data 9458b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points)); 9464b3e95d5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelData_Cuda_gen(ceed, num_input_fields, op_input_fields, qf_input_fields, num_output_fields, op_output_fields, 9474b3e95d5SJeremy L Thompson qf_output_fields, &data->max_P_1d, &Q_1d, &dim, &is_tensor, &use_3d_slices)); 9484b3e95d5SJeremy L Thompson if (dim == 0) dim = 1; 9494b3e95d5SJeremy L Thompson data->dim = dim; 9508b97b69aSJeremy L Thompson if (is_at_points) { 9518b97b69aSJeremy L Thompson CeedElemRestriction_Cuda *rstr_data; 9528b97b69aSJeremy L Thompson CeedElemRestriction rstr_points = NULL; 9534b3e95d5SJeremy L Thompson 9548b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL)); 9558b97b69aSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points)); 9568b97b69aSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetCompStride(rstr_points, &coords_comp_stride)); 9578b97b69aSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(rstr_points, &rstr_data)); 9588b97b69aSJeremy L Thompson data->points.indices = (CeedInt *)rstr_data->d_offsets; 9598b97b69aSJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points)); 9608b97b69aSJeremy L Thompson } 9618b97b69aSJeremy L Thompson if (is_at_points) use_3d_slices = false; 9628b97b69aSJeremy L Thompson if (Q_1d == 0) { 9638b97b69aSJeremy L Thompson if (is_at_points) Q_1d = max_num_points; 9648b97b69aSJeremy L Thompson else CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q_1d)); 9654b3e95d5SJeremy L Thompson } 9664b3e95d5SJeremy L Thompson data->Q_1d = Q_1d; 9674b3e95d5SJeremy L Thompson 9680b454692Sjeremylt // Check for restriction only identity operator 9694b3e95d5SJeremy L Thompson { 9704b3e95d5SJeremy L Thompson bool is_identity_qf; 9714b3e95d5SJeremy L Thompson 9722b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf)); 9730b454692Sjeremylt if (is_identity_qf) { 9749e201c85SYohann CeedEvalMode eval_mode_in, eval_mode_out; 975ca735530SJeremy L Thompson 9762b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in)); 9772b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out)); 9786574a04fSJeremy L Thompson CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND, 9796574a04fSJeremy L Thompson "Backend does not implement restriction only identity operators"); 9800b454692Sjeremylt } 9814b3e95d5SJeremy L Thompson } 9820b454692Sjeremylt 983f1a13f77SYohann Dudouit // Add atomicAdd function for old NVidia architectures 9844b3e95d5SJeremy L Thompson { 9854b3e95d5SJeremy L Thompson Ceed_Cuda *ceed_data; 9864b3e95d5SJeremy L Thompson struct cudaDeviceProp prop; 9874b3e95d5SJeremy L Thompson 9882b730f8bSJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_data)); 9892b730f8bSJeremy L Thompson CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id)); 99080a9ef05SNatalie Beams if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) { 9919c25dd66SJeremy L Thompson code << "// AtomicAdd fallback source\n"; 9929c25dd66SJeremy L Thompson code << "#include <ceed/jit-source/cuda/cuda-atomic-add-fallback.h>\n\n"; 993f1a13f77SYohann Dudouit } 9944b3e95d5SJeremy L Thompson } 995f1a13f77SYohann Dudouit 9969e201c85SYohann // Load basis source files 997dc007f05SJeremy L Thompson if (is_tensor) { 9989c25dd66SJeremy L Thompson code << "// Tensor basis source\n"; 9999c25dd66SJeremy L Thompson code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h>\n\n"; 1000dc007f05SJeremy L Thompson } else { 1001dc007f05SJeremy L Thompson code << "// Non-tensor basis source\n"; 1002dc007f05SJeremy L Thompson code << "#include <ceed/jit-source/cuda/cuda-shared-basis-nontensor-templates.h>\n\n"; 1003dc007f05SJeremy L Thompson } 1004dc007f05SJeremy L Thompson if (is_at_points) { 10058b97b69aSJeremy L Thompson code << "// AtPoints basis source\n"; 10068b97b69aSJeremy L Thompson code << "#include <ceed/jit-source/cuda/cuda-shared-basis-tensor-at-points-templates.h>\n\n"; 1007dc007f05SJeremy L Thompson } 10089c25dd66SJeremy L Thompson code << "// CodeGen operator source\n"; 10099c25dd66SJeremy L Thompson code << "#include <ceed/jit-source/cuda/cuda-gen-templates.h>\n\n"; 1010241a4b83SYohann 10114b3e95d5SJeremy L Thompson // Get QFunction name 10124b3e95d5SJeremy L Thompson std::string qfunction_name(qf_data->qfunction_name); 10134b3e95d5SJeremy L Thompson std::string operator_name; 10144b3e95d5SJeremy L Thompson 101509095acaSJeremy L Thompson operator_name = "CeedKernelCudaGenOperator_" + qfunction_name; 10164d537eeaSYohann 10179e201c85SYohann // Define CEED_Q_VLA 10189e201c85SYohann code << "\n#undef CEED_Q_VLA\n"; 1019dc007f05SJeremy L Thompson if (dim != 3 || is_at_points || use_3d_slices || !is_tensor) { 10209e201c85SYohann code << "#define CEED_Q_VLA 1\n\n"; 10219e201c85SYohann } else { 10229e201c85SYohann code << "#define CEED_Q_VLA " << Q_1d << "\n\n"; 10239e201c85SYohann } 10249e201c85SYohann 10254b3e95d5SJeremy L Thompson // Add user QFunction source 10264b3e95d5SJeremy L Thompson { 10279c25dd66SJeremy L Thompson const char *source_path; 10284b3e95d5SJeremy L Thompson 10299c25dd66SJeremy L Thompson CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path)); 10309c25dd66SJeremy L Thompson CeedCheck(source_path, ceed, CEED_ERROR_UNSUPPORTED, "/gpu/cuda/gen backend requires QFunction source code file"); 10319c25dd66SJeremy L Thompson 10329c25dd66SJeremy L Thompson code << "// User QFunction source\n"; 10339c25dd66SJeremy L Thompson code << "#include \"" << source_path << "\"\n\n"; 10344b3e95d5SJeremy L Thompson } 1035241a4b83SYohann 1036241a4b83SYohann // Setup 1037818e0025SJeremy L Thompson code << "\n// -----------------------------------------------------------------------------\n"; 10384b3e95d5SJeremy L Thompson code << "// Operator Kernel\n"; 10394b3e95d5SJeremy L Thompson code << "// \n"; 10404b3e95d5SJeremy L Thompson code << "// d_[in,out]_i: CeedVector device array\n"; 10414b3e95d5SJeremy L Thompson code << "// r_[in,out]_e_i: Element vector register\n"; 10424b3e95d5SJeremy L Thompson code << "// r_[in,out]_q_i: Quadrature space vector register\n"; 1043*9123fb08SJeremy L Thompson code << "// r_[in,out]_c_i: AtPoints Chebyshev coefficients register\n"; 10444b3e95d5SJeremy L Thompson code << "// r_[in,out]_s_i: Quadrature space slice vector register\n"; 10454b3e95d5SJeremy L Thompson code << "// \n"; 10464b3e95d5SJeremy L Thompson code << "// s_B_[in,out]_i: Interpolation matrix, shared memory\n"; 10474b3e95d5SJeremy L Thompson code << "// s_G_[in,out]_i: Gradient matrix, shared memory\n"; 10484b3e95d5SJeremy L Thompson code << "// -----------------------------------------------------------------------------\n"; 10494b3e95d5SJeremy L Thompson code << "extern \"C\" __global__ void " << operator_name 10508b97b69aSJeremy L Thompson << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar *W, Points_Cuda " 10518b97b69aSJeremy L Thompson "points) {\n"; 10524b3e95d5SJeremy L Thompson 10534b3e95d5SJeremy L Thompson // Scratch buffers 10549e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 10554b3e95d5SJeremy L Thompson CeedEvalMode eval_mode; 10564b3e95d5SJeremy L Thompson 10572b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 10589e201c85SYohann if (eval_mode != CEED_EVAL_WEIGHT) { // Skip CEED_EVAL_WEIGHT 10594b3e95d5SJeremy L Thompson code << " const CeedScalar *d_in_" << i << " = fields.inputs[" << i << "];\n"; 1060241a4b83SYohann } 1061241a4b83SYohann } 10629e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 10634b3e95d5SJeremy L Thompson code << " CeedScalar *d_out_" << i << " = fields.outputs[" << i << "];\n"; 1064241a4b83SYohann } 10651da99368SJeremy L Thompson 10669e201c85SYohann code << " const CeedInt dim = " << dim << ";\n"; 1067dc007f05SJeremy L Thompson code << " const CeedInt " << (is_tensor ? "Q_1d" : "Q") << " = " << Q_1d << ";\n"; 10688b97b69aSJeremy L Thompson if (is_at_points) { 10698b97b69aSJeremy L Thompson code << " const CeedInt max_num_points = " << max_num_points << ";\n"; 10708b97b69aSJeremy L Thompson code << " const CeedInt coords_comp_stride = " << coords_comp_stride << ";\n"; 10718b97b69aSJeremy L Thompson } 10721da99368SJeremy L Thompson 10734b3e95d5SJeremy L Thompson // Shared data 1074241a4b83SYohann code << " extern __shared__ CeedScalar slice[];\n"; 10759e201c85SYohann code << " SharedData_Cuda data;\n"; 10769e201c85SYohann code << " data.t_id_x = threadIdx.x;\n"; 10779e201c85SYohann code << " data.t_id_y = threadIdx.y;\n"; 10789e201c85SYohann code << " data.t_id_z = threadIdx.z;\n"; 10799e201c85SYohann code << " data.t_id = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n"; 1080dc007f05SJeremy L Thompson code << " data.slice = slice + data.t_id_z*T_1D" << ((!is_tensor || dim == 1) ? "" : "*T_1D") << ";\n"; 1081920dcdc4Sjeremylt 1082ac421f39SYohann // Initialize constants, and matrices B and G 10834b3e95d5SJeremy L Thompson code << "\n // Input field constants and basis data\n"; 10849e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 1085dc007f05SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_input_fields[i], qf_input_fields[i], Q_1d, true, is_tensor, 1086dc007f05SJeremy L Thompson is_at_points, use_3d_slices)); 1087920dcdc4Sjeremylt } 10884b3e95d5SJeremy L Thompson code << "\n // Output field constants and basis data\n"; 10899e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 1090dc007f05SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelFieldData_Cuda_gen(code, data, i, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_tensor, 1091dc007f05SJeremy L Thompson is_at_points, use_3d_slices)); 10924b3e95d5SJeremy L Thompson } 1093920dcdc4Sjeremylt 10944b3e95d5SJeremy L Thompson // Loop over all elements 10954b3e95d5SJeremy L Thompson code << "\n // Element loop\n"; 1096ac421f39SYohann code << " __syncthreads();\n"; 10979e201c85SYohann code << " for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n"; 10984b3e95d5SJeremy L Thompson 1099e93651e5SJeremy L Thompson // -- Compute minimum buffer space needed 11008b97b69aSJeremy L Thompson CeedInt max_rstr_buffer_size = 1; 1101e93651e5SJeremy L Thompson 11029e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 1103e93651e5SJeremy L Thompson CeedInt num_comp, elem_size; 1104e93651e5SJeremy L Thompson CeedElemRestriction elem_rstr; 1105e93651e5SJeremy L Thompson 1106e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 1107e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 1108e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1109dc007f05SJeremy L Thompson max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1)); 1110681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 1111e93651e5SJeremy L Thompson } 1112e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 1113e93651e5SJeremy L Thompson CeedInt num_comp, elem_size; 1114e93651e5SJeremy L Thompson CeedElemRestriction elem_rstr; 1115e93651e5SJeremy L Thompson 1116e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 1117e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 1118e93651e5SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1119dc007f05SJeremy L Thompson max_rstr_buffer_size = CeedIntMax(max_rstr_buffer_size, num_comp * (is_tensor && (dim >= 3) ? elem_size : 1)); 1120681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 1121e93651e5SJeremy L Thompson } 1122e93651e5SJeremy L Thompson code << " // Scratch restriction buffer space\n"; 1123e93651e5SJeremy L Thompson code << " CeedScalar r_e_scratch[" << max_rstr_buffer_size << "];\n"; 1124e93651e5SJeremy L Thompson 1125e93651e5SJeremy L Thompson // -- Determine best input field processing order 1126e93651e5SJeremy L Thompson CeedInt field_rstr_in_buffer[CEED_FIELD_MAX], input_field_order[CEED_FIELD_MAX]; 1127e93651e5SJeremy L Thompson 1128e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 1129e93651e5SJeremy L Thompson field_rstr_in_buffer[i] = -1; 1130e93651e5SJeremy L Thompson input_field_order[i] = -1; 1131e93651e5SJeremy L Thompson } 1132e93651e5SJeremy L Thompson { 1133e93651e5SJeremy L Thompson bool is_ordered[CEED_FIELD_MAX]; 1134e93651e5SJeremy L Thompson CeedInt curr_index = 0; 1135e93651e5SJeremy L Thompson 1136e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) is_ordered[i] = false; 1137e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 1138e93651e5SJeremy L Thompson CeedVector vec_i; 1139e93651e5SJeremy L Thompson CeedElemRestriction rstr_i; 1140e93651e5SJeremy L Thompson 1141e93651e5SJeremy L Thompson if (is_ordered[i]) continue; 1142e93651e5SJeremy L Thompson field_rstr_in_buffer[i] = i; 1143e93651e5SJeremy L Thompson is_ordered[i] = true; 1144e93651e5SJeremy L Thompson input_field_order[curr_index] = i; 1145e93651e5SJeremy L Thompson curr_index++; 1146034f99fdSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i)); 1147e93651e5SJeremy L Thompson if (vec_i == CEED_VECTOR_NONE) continue; // CEED_EVAL_WEIGHT 1148e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i)); 1149e93651e5SJeremy L Thompson for (CeedInt j = i + 1; j < num_input_fields; j++) { 1150e93651e5SJeremy L Thompson CeedVector vec_j; 1151e93651e5SJeremy L Thompson CeedElemRestriction rstr_j; 1152e93651e5SJeremy L Thompson 1153e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[j], &vec_j)); 1154e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[j], &rstr_j)); 1155e93651e5SJeremy L Thompson if (rstr_i == rstr_j && vec_i == vec_j) { 1156e93651e5SJeremy L Thompson field_rstr_in_buffer[j] = i; 1157e93651e5SJeremy L Thompson is_ordered[j] = true; 1158e93651e5SJeremy L Thompson input_field_order[curr_index] = j; 1159e93651e5SJeremy L Thompson curr_index++; 1160e93651e5SJeremy L Thompson } 1161681d0ea7SJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec_j)); 1162681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j)); 1163e93651e5SJeremy L Thompson } 1164681d0ea7SJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec_i)); 1165681d0ea7SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i)); 1166e93651e5SJeremy L Thompson } 1167e93651e5SJeremy L Thompson } 1168e93651e5SJeremy L Thompson 1169e93651e5SJeremy L Thompson // -- Input restriction and basis 1170e93651e5SJeremy L Thompson code << "\n // -- Input field restrictions and basis actions\n"; 1171e93651e5SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 1172e93651e5SJeremy L Thompson CeedInt f = input_field_order[i]; 1173e93651e5SJeremy L Thompson 1174e93651e5SJeremy L Thompson code << " // ---- Input field " << f << "\n"; 1175920dcdc4Sjeremylt 11764b3e95d5SJeremy L Thompson // ---- Restriction 1177e93651e5SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, f, dim, field_rstr_in_buffer, op_input_fields[f], qf_input_fields[f], 1178dc007f05SJeremy L Thompson Q_1d, true, is_tensor, is_at_points, use_3d_slices)); 11799a2291e3SJeremy L Thompson 11804b3e95d5SJeremy L Thompson // ---- Basis action 1181dc007f05SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, f, dim, op_input_fields[f], qf_input_fields[f], Q_1d, true, is_tensor, 1182dc007f05SJeremy L Thompson is_at_points, use_3d_slices)); 1183920dcdc4Sjeremylt } 1184920dcdc4Sjeremylt 11854b3e95d5SJeremy L Thompson // -- Q function 11868b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelQFunction_Cuda_gen(code, data, dim, max_num_points, num_input_fields, op_input_fields, qf_input_fields, 1187dc007f05SJeremy L Thompson num_output_fields, op_output_fields, qf_output_fields, qfunction_name, Q_1d, is_tensor, 1188dc007f05SJeremy L Thompson is_at_points, use_3d_slices)); 1189ca735530SJeremy L Thompson 11904b3e95d5SJeremy L Thompson // -- Output basis and restriction 11914b3e95d5SJeremy L Thompson code << "\n // -- Output field basis action and restrictions\n"; 11929e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 11934b3e95d5SJeremy L Thompson code << " // ---- Output field " << i << "\n"; 1194ca735530SJeremy L Thompson 11954b3e95d5SJeremy L Thompson // ---- Basis action 1196dc007f05SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelBasis_Cuda_gen(code, data, i, dim, op_output_fields[i], qf_output_fields[i], Q_1d, false, is_tensor, 1197dc007f05SJeremy L Thompson is_at_points, use_3d_slices)); 11989a2291e3SJeremy L Thompson 11994b3e95d5SJeremy L Thompson // ---- Restriction 12008b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernelRestriction_Cuda_gen(code, data, i, dim, NULL, op_output_fields[i], qf_output_fields[i], Q_1d, false, 1201dc007f05SJeremy L Thompson is_tensor, is_at_points, use_3d_slices)); 1202ac421f39SYohann } 1203241a4b83SYohann 12044b3e95d5SJeremy L Thompson // Close loop and function 1205241a4b83SYohann code << " }\n"; 1206818e0025SJeremy L Thompson code << "}\n"; 1207818e0025SJeremy L Thompson code << "// -----------------------------------------------------------------------------\n\n"; 1208241a4b83SYohann 12093a2968d6SJeremy L Thompson // Compile 1210eb7e6cafSJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, code.str().c_str(), &data->module, 1, "T_1D", CeedIntMax(Q_1d, data->max_P_1d))); 1211eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op)); 12122b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorSetSetupDone(op)); 12139bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 1214c11e12f4SJeremy L Thompson CeedCallBackend(CeedQFunctionDestroy(&qf)); 1215e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1216241a4b83SYohann } 12172a86cc9dSSebastian Grimberg 1218ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 1219