1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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
1605c335cbSnbeams #include "../hip/ceed-hip-common.h"
170d0321e0SJeremy L Thompson #include "../hip/ceed-hip-compile.h"
182b730f8bSJeremy L Thompson #include "ceed-hip-ref.h"
190d0321e0SJeremy L Thompson
200d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
210d0321e0SJeremy L Thompson // Build QFunction kernel
220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
CeedQFunctionBuildKernel_Hip_ref(CeedQFunction qf)23eb7e6cafSJeremy L Thompson extern "C" int CeedQFunctionBuildKernel_Hip_ref(CeedQFunction qf) {
240d0321e0SJeremy L Thompson using std::ostringstream;
250d0321e0SJeremy L Thompson using std::string;
26b7453713SJeremy L Thompson
27437930d1SJeremy L Thompson Ceed ceed;
2886e1ed65Snbeams Ceed_Hip *ceed_Hip;
29b7453713SJeremy L Thompson CeedInt num_input_fields, num_output_fields, size;
30b7453713SJeremy L Thompson CeedQFunctionField *input_fields, *output_fields;
310d0321e0SJeremy L Thompson CeedQFunction_Hip *data;
32b7453713SJeremy L Thompson
339bc66399SJeremy L Thompson // QFunction is built
349bc66399SJeremy L Thompson CeedCallBackend(CeedQFunctionGetData(qf, (void **)&data));
359bc66399SJeremy L Thompson if (data->QFunction) return CEED_ERROR_SUCCESS;
369bc66399SJeremy L Thompson
376e536b99SJeremy L Thompson CeedCallBackend(CeedQFunctionGetCeed(qf, &ceed));
38b7453713SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &ceed_Hip));
390d0321e0SJeremy L Thompson
400d0321e0SJeremy L Thompson // QFunction kernel generation
412b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionGetFields(qf, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
420d0321e0SJeremy L Thompson
430d0321e0SJeremy L Thompson // Build strings for final kernel
44437930d1SJeremy L Thompson string qfunction_name(data->qfunction_name);
45204bfdd7SJeremy L Thompson string kernel_name = "CeedKernelHipRefQFunction_" + qfunction_name;
460d0321e0SJeremy L Thompson ostringstream code;
470d0321e0SJeremy L Thompson
489c25dd66SJeremy L Thompson // Definitions
499c25dd66SJeremy L Thompson code << "// QFunction source\n";
509c25dd66SJeremy L Thompson code << "#include <ceed/jit-source/hip/hip-ref-qfunction.h>\n\n";
519c25dd66SJeremy L Thompson {
529c25dd66SJeremy L Thompson const char *source_path;
539c25dd66SJeremy L Thompson
549c25dd66SJeremy L Thompson CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
559c25dd66SJeremy L Thompson CeedCheck(source_path, ceed, CEED_ERROR_BACKEND, "No QFunction source or hipFunction_t provided.");
569c25dd66SJeremy L Thompson
579c25dd66SJeremy L Thompson code << "// User QFunction source\n";
589c25dd66SJeremy L Thompson code << "#include \"" << source_path << "\"\n\n";
599c25dd66SJeremy L Thompson }
6086e1ed65Snbeams code << "extern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
6186e1ed65Snbeams code << "__global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Hip fields) {\n";
620d0321e0SJeremy L Thompson
630d0321e0SJeremy L Thompson // Inputs
6446dc0734SJeremy L Thompson code << " // Input fields\n";
65437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) {
662b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetSize(input_fields[i], &size));
6746dc0734SJeremy L Thompson code << " const CeedInt size_input_" << i << " = " << size << ";\n";
6846dc0734SJeremy L Thompson code << " CeedScalar input_" << i << "[size_input_" << i << "];\n";
6946dc0734SJeremy L Thompson }
709b443e3bSJeremy L Thompson code << " const CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
7146dc0734SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) {
7246dc0734SJeremy L Thompson code << " inputs[" << i << "] = input_" << i << ";\n";
730d0321e0SJeremy L Thompson }
74437930d1SJeremy L Thompson code << "\n";
750d0321e0SJeremy L Thompson
760d0321e0SJeremy L Thompson // Outputs
7746dc0734SJeremy L Thompson code << " // Output fields\n";
78437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) {
792b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetSize(output_fields[i], &size));
8046dc0734SJeremy L Thompson code << " const CeedInt size_output_" << i << " = " << size << ";\n";
8146dc0734SJeremy L Thompson code << " CeedScalar output_" << i << "[size_output_" << i << "];\n";
820d0321e0SJeremy L Thompson }
839b443e3bSJeremy L Thompson code << " CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
84437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) {
8546dc0734SJeremy L Thompson code << " outputs[" << i << "] = output_" << i << ";\n";
860d0321e0SJeremy L Thompson }
87437930d1SJeremy L Thompson code << "\n";
880d0321e0SJeremy L Thompson
890d0321e0SJeremy L Thompson // Loop over quadrature points
9046dc0734SJeremy L Thompson code << " // Loop over quadrature points\n";
910d0321e0SJeremy L Thompson code << " for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n";
920d0321e0SJeremy L Thompson
930d0321e0SJeremy L Thompson // Load inputs
9446dc0734SJeremy L Thompson code << " // -- Load inputs\n";
95437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) {
9646dc0734SJeremy L Thompson code << " readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n";
970d0321e0SJeremy L Thompson }
9846dc0734SJeremy L Thompson code << "\n";
9946dc0734SJeremy L Thompson
1000d0321e0SJeremy L Thompson // QFunction
10146dc0734SJeremy L Thompson code << " // -- Call QFunction\n";
10246dc0734SJeremy L Thompson code << " " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n";
1030d0321e0SJeremy L Thompson
1040d0321e0SJeremy L Thompson // Write outputs
10546dc0734SJeremy L Thompson code << " // -- Write outputs\n";
106437930d1SJeremy L Thompson for (CeedInt i = 0; i < num_output_fields; i++) {
10746dc0734SJeremy L Thompson code << " writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n";
1080d0321e0SJeremy L Thompson }
1090d0321e0SJeremy L Thompson code << " }\n";
1100d0321e0SJeremy L Thompson code << "}\n";
1110d0321e0SJeremy L Thompson
1120d0321e0SJeremy L Thompson // Compile kernel
113eb7e6cafSJeremy L Thompson CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 1, "BLOCK_SIZE", ceed_Hip->opt_block_size));
114eb7e6cafSJeremy L Thompson CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, kernel_name.c_str(), &data->QFunction));
1159bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed));
1160d0321e0SJeremy L Thompson return CEED_ERROR_SUCCESS;
1170d0321e0SJeremy L Thompson }
1182a86cc9dSSebastian Grimberg
1190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
120