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. 30d0321e0SJeremy L Thompson // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 50d0321e0SJeremy L Thompson // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 70d0321e0SJeremy L Thompson 849aac155SJeremy L Thompson #include <ceed.h> 90d0321e0SJeremy L Thompson #include <ceed/backend.h> 10437930d1SJeremy L Thompson #include <ceed/jit-tools.h> 112b730f8bSJeremy L Thompson #include <string.h> 122b730f8bSJeremy L Thompson 130d0321e0SJeremy L Thompson #include <iostream> 140d0321e0SJeremy L Thompson #include <sstream> 152b730f8bSJeremy L Thompson 160d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h" 172b730f8bSJeremy L Thompson #include "ceed-cuda-ref.h" 180d0321e0SJeremy L Thompson 190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 200d0321e0SJeremy L Thompson // Build QFunction kernel 210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 22eb7e6cafSJeremy L Thompson extern "C" int CeedQFunctionBuildKernel_Cuda_ref(CeedQFunction qf) { 230d0321e0SJeremy L Thompson using std::ostringstream; 240d0321e0SJeremy L Thompson using std::string; 25ca735530SJeremy L Thompson 260d0321e0SJeremy L Thompson Ceed ceed; 27ca735530SJeremy L Thompson CeedInt num_input_fields, num_output_fields, size; 28ca735530SJeremy L Thompson CeedQFunctionField *input_fields, *output_fields; 290d0321e0SJeremy L Thompson CeedQFunction_Cuda *data; 30ca735530SJeremy L Thompson 310d0321e0SJeremy L Thompson // QFunction is built 32*9bc66399SJeremy L Thompson CeedCallBackend(CeedQFunctionGetData(qf, (void **)&data)); 332b730f8bSJeremy L Thompson if (data->QFunction) return CEED_ERROR_SUCCESS; 34437930d1SJeremy L Thompson 35*9bc66399SJeremy L Thompson CeedCallBackend(CeedQFunctionGetCeed(qf, &ceed)); 36*9bc66399SJeremy L Thompson 370d0321e0SJeremy L Thompson // QFunction kernel generation 382b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionGetFields(qf, &num_input_fields, &input_fields, &num_output_fields, &output_fields)); 390d0321e0SJeremy L Thompson 400d0321e0SJeremy L Thompson // Build strings for final kernel 41437930d1SJeremy L Thompson string qfunction_name(data->qfunction_name); 42204bfdd7SJeremy L Thompson string kernel_name = "CeedKernelCudaRefQFunction_" + qfunction_name; 430d0321e0SJeremy L Thompson ostringstream code; 440d0321e0SJeremy L Thompson 459c25dd66SJeremy L Thompson // Definitions 469c25dd66SJeremy L Thompson code << "// QFunction source\n"; 479c25dd66SJeremy L Thompson code << "#include <ceed/jit-source/cuda/cuda-ref-qfunction.h>\n\n"; 489c25dd66SJeremy L Thompson { 499c25dd66SJeremy L Thompson const char *source_path; 509c25dd66SJeremy L Thompson 519c25dd66SJeremy L Thompson CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path)); 529c25dd66SJeremy L Thompson CeedCheck(source_path, ceed, CEED_ERROR_BACKEND, "No QFunction source or CUfunction provided."); 539c25dd66SJeremy L Thompson 549c25dd66SJeremy L Thompson code << "// User QFunction source\n"; 559c25dd66SJeremy L Thompson code << "#include \"" << source_path << "\"\n\n"; 569c25dd66SJeremy L Thompson } 57437930d1SJeremy L Thompson code << "extern \"C\" __global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Cuda fields) {\n"; 580d0321e0SJeremy L Thompson 590d0321e0SJeremy L Thompson // Inputs 6046dc0734SJeremy L Thompson code << " // Input fields\n"; 61437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 622b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetSize(input_fields[i], &size)); 6346dc0734SJeremy L Thompson code << " const CeedInt size_input_" << i << " = " << size << ";\n"; 6446dc0734SJeremy L Thompson code << " CeedScalar input_" << i << "[size_input_" << i << "];\n"; 6546dc0734SJeremy L Thompson } 669b443e3bSJeremy L Thompson code << " const CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n"; 6746dc0734SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 6846dc0734SJeremy L Thompson code << " inputs[" << i << "] = input_" << i << ";\n"; 690d0321e0SJeremy L Thompson } 70437930d1SJeremy L Thompson code << "\n"; 710d0321e0SJeremy L Thompson 720d0321e0SJeremy L Thompson // Outputs 7346dc0734SJeremy L Thompson code << " // Output fields\n"; 74437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 752b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetSize(output_fields[i], &size)); 7646dc0734SJeremy L Thompson code << " const CeedInt size_output_" << i << " = " << size << ";\n"; 7746dc0734SJeremy L Thompson code << " CeedScalar output_" << i << "[size_output_" << i << "];\n"; 780d0321e0SJeremy L Thompson } 799b443e3bSJeremy L Thompson code << " CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n"; 80437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 8146dc0734SJeremy L Thompson code << " outputs[" << i << "] = output_" << i << ";\n"; 820d0321e0SJeremy L Thompson } 83437930d1SJeremy L Thompson code << "\n"; 840d0321e0SJeremy L Thompson 850d0321e0SJeremy L Thompson // Loop over quadrature points 8646dc0734SJeremy L Thompson code << " // Loop over quadrature points\n"; 870d0321e0SJeremy L Thompson code << " for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n"; 880d0321e0SJeremy L Thompson 890d0321e0SJeremy L Thompson // Load inputs 9046dc0734SJeremy L Thompson code << " // -- Load inputs\n"; 91437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 9246dc0734SJeremy L Thompson code << " readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n"; 930d0321e0SJeremy L Thompson } 9446dc0734SJeremy L Thompson code << "\n"; 9546dc0734SJeremy L Thompson 960d0321e0SJeremy L Thompson // QFunction 9746dc0734SJeremy L Thompson code << " // -- Call QFunction\n"; 9846dc0734SJeremy L Thompson code << " " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n"; 990d0321e0SJeremy L Thompson 1000d0321e0SJeremy L Thompson // Write outputs 10146dc0734SJeremy L Thompson code << " // -- Write outputs\n"; 102437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) { 10346dc0734SJeremy L Thompson code << " writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n"; 1040d0321e0SJeremy L Thompson } 1050d0321e0SJeremy L Thompson code << " }\n"; 1060d0321e0SJeremy L Thompson code << "}\n"; 1070d0321e0SJeremy L Thompson 1080d0321e0SJeremy L Thompson // View kernel for debugging 10923d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated QFunction Kernels:\n"); 1100d0321e0SJeremy L Thompson CeedDebug(ceed, code.str().c_str()); 1110d0321e0SJeremy L Thompson 1120d0321e0SJeremy L Thompson // Compile kernel 113eb7e6cafSJeremy L Thompson CeedCallBackend(CeedCompile_Cuda(ceed, code.str().c_str(), &data->module, 0)); 114eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, kernel_name.c_str(), &data->QFunction)); 115*9bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 1160d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS; 1170d0321e0SJeremy L Thompson } 1182a86cc9dSSebastian Grimberg 1190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 120