xref: /libCEED/backends/hip-ref/ceed-hip-ref-qfunction-load.cpp (revision 23d4529ec201513b6e781993a7a761deb8126e79)
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 
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 //------------------------------------------------------------------------------
23eb7e6cafSJeremy L Thompson extern "C" int CeedQFunctionBuildKernel_Hip_ref(CeedQFunction qf) {
240d0321e0SJeremy L Thompson   using std::ostringstream;
250d0321e0SJeremy L Thompson   using std::string;
26437930d1SJeremy L Thompson   Ceed ceed;
27437930d1SJeremy L Thompson   CeedQFunctionGetCeed(qf, &ceed);
2886e1ed65Snbeams   Ceed_Hip *ceed_Hip;
292b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &ceed_Hip));
300d0321e0SJeremy L Thompson   CeedQFunction_Hip *data;
312b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, (void **)&data));
32437930d1SJeremy L Thompson 
330d0321e0SJeremy L Thompson   // QFunction is built
342b730f8bSJeremy L Thompson   if (data->QFunction) return CEED_ERROR_SUCCESS;
350d0321e0SJeremy L Thompson 
366574a04fSJeremy L Thompson   CeedCheck(data->qfunction_source, ceed, CEED_ERROR_BACKEND, "No QFunction source or hipFunction_t provided.");
37437930d1SJeremy L Thompson 
380d0321e0SJeremy L Thompson   // QFunction kernel generation
39437930d1SJeremy L Thompson   CeedInt             num_input_fields, num_output_fields, size;
40437930d1SJeremy L Thompson   CeedQFunctionField *input_fields, *output_fields;
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   char *read_write_kernel_path, *read_write_kernel_source;
452b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-ref-qfunction.h", &read_write_kernel_path));
46*23d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading QFunction Read/Write Kernel Source -----\n");
472b730f8bSJeremy L Thompson   CeedCallBackend(CeedLoadSourceToBuffer(ceed, read_write_kernel_path, &read_write_kernel_source));
48*23d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading QFunction Read/Write Kernel Source Complete! -----\n");
49437930d1SJeremy L Thompson   string        qfunction_source(data->qfunction_source);
50437930d1SJeremy L Thompson   string        qfunction_name(data->qfunction_name);
51437930d1SJeremy L Thompson   string        read_write(read_write_kernel_source);
52204bfdd7SJeremy L Thompson   string        kernel_name = "CeedKernelHipRefQFunction_" + qfunction_name;
530d0321e0SJeremy L Thompson   ostringstream code;
540d0321e0SJeremy L Thompson 
550d0321e0SJeremy L Thompson   // Defintions
56437930d1SJeremy L Thompson   code << read_write;
57437930d1SJeremy L Thompson   code << qfunction_source;
5846dc0734SJeremy L Thompson   code << "\n";
5986e1ed65Snbeams   code << "extern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
6086e1ed65Snbeams   code << "__global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Hip fields) {\n";
610d0321e0SJeremy L Thompson 
620d0321e0SJeremy L Thompson   // Inputs
6346dc0734SJeremy L Thompson   code << "  // Input fields\n";
64437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
652b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(input_fields[i], &size));
6646dc0734SJeremy L Thompson     code << "  const CeedInt size_input_" << i << " = " << size << ";\n";
6746dc0734SJeremy L Thompson     code << "  CeedScalar input_" << i << "[size_input_" << i << "];\n";
6846dc0734SJeremy L Thompson   }
6946dc0734SJeremy L Thompson   code << "  const CeedScalar* inputs[" << num_input_fields << "];\n";
7046dc0734SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
7146dc0734SJeremy L Thompson     code << "  inputs[" << i << "] = input_" << i << ";\n";
720d0321e0SJeremy L Thompson   }
73437930d1SJeremy L Thompson   code << "\n";
740d0321e0SJeremy L Thompson 
750d0321e0SJeremy L Thompson   // Outputs
7646dc0734SJeremy L Thompson   code << "  // Output fields\n";
77437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
782b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetSize(output_fields[i], &size));
7946dc0734SJeremy L Thompson     code << "  const CeedInt size_output_" << i << " = " << size << ";\n";
8046dc0734SJeremy L Thompson     code << "  CeedScalar output_" << i << "[size_output_" << i << "];\n";
810d0321e0SJeremy L Thompson   }
8246dc0734SJeremy L Thompson   code << "  CeedScalar* outputs[" << num_output_fields << "];\n";
83437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
8446dc0734SJeremy L Thompson     code << "  outputs[" << i << "] = output_" << i << ";\n";
850d0321e0SJeremy L Thompson   }
86437930d1SJeremy L Thompson   code << "\n";
870d0321e0SJeremy L Thompson 
880d0321e0SJeremy L Thompson   // Loop over quadrature points
8946dc0734SJeremy L Thompson   code << "  // Loop over quadrature points\n";
900d0321e0SJeremy L Thompson   code << "  for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n";
910d0321e0SJeremy L Thompson 
920d0321e0SJeremy L Thompson   // Load inputs
9346dc0734SJeremy L Thompson   code << "    // -- Load inputs\n";
94437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
9546dc0734SJeremy L Thompson     code << "    readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n";
960d0321e0SJeremy L Thompson   }
9746dc0734SJeremy L Thompson   code << "\n";
9846dc0734SJeremy L Thompson 
990d0321e0SJeremy L Thompson   // QFunction
10046dc0734SJeremy L Thompson   code << "    // -- Call QFunction\n";
10146dc0734SJeremy L Thompson   code << "    " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n";
1020d0321e0SJeremy L Thompson 
1030d0321e0SJeremy L Thompson   // Write outputs
10446dc0734SJeremy L Thompson   code << "    // -- Write outputs\n";
105437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
10646dc0734SJeremy L Thompson     code << "    writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n";
1070d0321e0SJeremy L Thompson   }
1080d0321e0SJeremy L Thompson   code << "  }\n";
1090d0321e0SJeremy L Thompson   code << "}\n";
1100d0321e0SJeremy L Thompson 
1110d0321e0SJeremy L Thompson   // View kernel for debugging
112*23d4529eSJeremy L Thompson   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated QFunction Kernels:\n");
1130d0321e0SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
1140d0321e0SJeremy L Thompson 
1150d0321e0SJeremy L Thompson   // Compile kernel
116eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 1, "BLOCK_SIZE", ceed_Hip->opt_block_size));
117eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, kernel_name.c_str(), &data->QFunction));
1180d0321e0SJeremy L Thompson 
1190d0321e0SJeremy L Thompson   // Cleanup
1202b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&data->qfunction_source));
1212b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&read_write_kernel_path));
1222b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&read_write_kernel_source));
123437930d1SJeremy L Thompson 
1240d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1250d0321e0SJeremy L Thompson }
1262a86cc9dSSebastian Grimberg 
1270d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
128