xref: /libCEED/backends/hip-ref/ceed-hip-ref-qfunction-load.cpp (revision 2b730f8b5a9c809740a0b3b302db43a719c636b1)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
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 
80d0321e0SJeremy L Thompson #include <ceed/backend.h>
9*2b730f8bSJeremy L Thompson #include <ceed/ceed.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
11*2b730f8bSJeremy L Thompson #include <string.h>
12*2b730f8bSJeremy L Thompson 
130d0321e0SJeremy L Thompson #include <iostream>
140d0321e0SJeremy L Thompson #include <sstream>
15*2b730f8bSJeremy L Thompson 
160d0321e0SJeremy L Thompson #include "../hip/ceed-hip-compile.h"
17*2b730f8bSJeremy L Thompson #include "ceed-hip-ref.h"
180d0321e0SJeremy L Thompson 
190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
200d0321e0SJeremy L Thompson // Build QFunction kernel
210d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
220d0321e0SJeremy L Thompson extern "C" int CeedHipBuildQFunction(CeedQFunction qf) {
230d0321e0SJeremy L Thompson   using std::ostringstream;
240d0321e0SJeremy L Thompson   using std::string;
25437930d1SJeremy L Thompson   Ceed ceed;
26437930d1SJeremy L Thompson   CeedQFunctionGetCeed(qf, &ceed);
2786e1ed65Snbeams   Ceed_Hip *ceed_Hip;
28*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Hip));
290d0321e0SJeremy L Thompson   CeedQFunction_Hip *data;
30*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, (void **)&data));
31437930d1SJeremy L Thompson 
320d0321e0SJeremy L Thompson   // QFunction is built
33*2b730f8bSJeremy L Thompson   if (data->QFunction) return CEED_ERROR_SUCCESS;
340d0321e0SJeremy L Thompson 
35*2b730f8bSJeremy L Thompson   if (!data->qfunction_source) {
36437930d1SJeremy L Thompson     // LCOV_EXCL_START
37*2b730f8bSJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND, "No QFunction source or hipFunction_t provided.");
38437930d1SJeremy L Thompson     // LCOV_EXCL_STOP
39*2b730f8bSJeremy L Thompson   }
40437930d1SJeremy L Thompson 
410d0321e0SJeremy L Thompson   // QFunction kernel generation
42437930d1SJeremy L Thompson   CeedInt             num_input_fields, num_output_fields, size;
43437930d1SJeremy L Thompson   CeedQFunctionField *input_fields, *output_fields;
44*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
450d0321e0SJeremy L Thompson 
460d0321e0SJeremy L Thompson   // Build strings for final kernel
47437930d1SJeremy L Thompson   char *read_write_kernel_path, *read_write_kernel_source;
48*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-ref-qfunction.h", &read_write_kernel_path));
4946dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading QFunction Read/Write Kernel Source -----\n");
50*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, read_write_kernel_path, &read_write_kernel_source));
5146dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading QFunction Read/Write Kernel Source Complete! -----\n");
52437930d1SJeremy L Thompson   string        qfunction_source(data->qfunction_source);
53437930d1SJeremy L Thompson   string        qfunction_name(data->qfunction_name);
54437930d1SJeremy L Thompson   string        read_write(read_write_kernel_source);
55204bfdd7SJeremy L Thompson   string        kernel_name = "CeedKernelHipRefQFunction_" + qfunction_name;
560d0321e0SJeremy L Thompson   ostringstream code;
570d0321e0SJeremy L Thompson 
580d0321e0SJeremy L Thompson   // Defintions
59437930d1SJeremy L Thompson   code << read_write;
60437930d1SJeremy L Thompson   code << qfunction_source;
6146dc0734SJeremy L Thompson   code << "\n";
6286e1ed65Snbeams   code << "extern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
6386e1ed65Snbeams   code << "__global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Hip fields) {\n";
640d0321e0SJeremy L Thompson 
650d0321e0SJeremy L Thompson   // Inputs
6646dc0734SJeremy L Thompson   code << "  // Input fields\n";
67437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
68*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(input_fields[i], &size));
6946dc0734SJeremy L Thompson     code << "  const CeedInt size_input_" << i << " = " << size << ";\n";
7046dc0734SJeremy L Thompson     code << "  CeedScalar input_" << i << "[size_input_" << i << "];\n";
7146dc0734SJeremy L Thompson   }
7246dc0734SJeremy L Thompson   code << "  const CeedScalar* inputs[" << num_input_fields << "];\n";
7346dc0734SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
7446dc0734SJeremy L Thompson     code << "  inputs[" << i << "] = input_" << i << ";\n";
750d0321e0SJeremy L Thompson   }
76437930d1SJeremy L Thompson   code << "\n";
770d0321e0SJeremy L Thompson 
780d0321e0SJeremy L Thompson   // Outputs
7946dc0734SJeremy L Thompson   code << "  // Output fields\n";
80437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
81*2b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(output_fields[i], &size));
8246dc0734SJeremy L Thompson     code << "  const CeedInt size_output_" << i << " = " << size << ";\n";
8346dc0734SJeremy L Thompson     code << "  CeedScalar output_" << i << "[size_output_" << i << "];\n";
840d0321e0SJeremy L Thompson   }
8546dc0734SJeremy L Thompson   code << "  CeedScalar* outputs[" << num_output_fields << "];\n";
86437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
8746dc0734SJeremy L Thompson     code << "  outputs[" << i << "] = output_" << i << ";\n";
880d0321e0SJeremy L Thompson   }
89437930d1SJeremy L Thompson   code << "\n";
900d0321e0SJeremy L Thompson 
910d0321e0SJeremy L Thompson   // Loop over quadrature points
9246dc0734SJeremy L Thompson   code << "  // Loop over quadrature points\n";
930d0321e0SJeremy L Thompson   code << "  for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n";
940d0321e0SJeremy L Thompson 
950d0321e0SJeremy L Thompson   // Load inputs
9646dc0734SJeremy L Thompson   code << "    // -- Load inputs\n";
97437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
9846dc0734SJeremy L Thompson     code << "    readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n";
990d0321e0SJeremy L Thompson   }
10046dc0734SJeremy L Thompson   code << "\n";
10146dc0734SJeremy L Thompson 
1020d0321e0SJeremy L Thompson   // QFunction
10346dc0734SJeremy L Thompson   code << "    // -- Call QFunction\n";
10446dc0734SJeremy L Thompson   code << "    " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n";
1050d0321e0SJeremy L Thompson 
1060d0321e0SJeremy L Thompson   // Write outputs
10746dc0734SJeremy L Thompson   code << "    // -- Write outputs\n";
108437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
10946dc0734SJeremy L Thompson     code << "    writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n";
1100d0321e0SJeremy L Thompson   }
1110d0321e0SJeremy L Thompson   code << "  }\n";
1120d0321e0SJeremy L Thompson   code << "}\n";
1130d0321e0SJeremy L Thompson 
1140d0321e0SJeremy L Thompson   // View kernel for debugging
11546dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "Generated QFunction Kernels:\n");
1160d0321e0SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
1170d0321e0SJeremy L Thompson 
1180d0321e0SJeremy L Thompson   // Compile kernel
119*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedCompileHip(ceed, code.str().c_str(), &data->module, 1, "BLOCK_SIZE", ceed_Hip->opt_block_size));
120*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetKernelHip(ceed, data->module, kernel_name.c_str(), &data->QFunction));
1210d0321e0SJeremy L Thompson 
1220d0321e0SJeremy L Thompson   // Cleanup
123*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data->qfunction_source));
124*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&read_write_kernel_path));
125*2b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&read_write_kernel_source));
126437930d1SJeremy L Thompson 
1270d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1280d0321e0SJeremy L Thompson }
1290d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
130