xref: /libCEED/backends/hip-ref/ceed-hip-ref-qfunction-load.cpp (revision 204bfdd766dd3467cc980d6ab2360b045ce7fa2c)
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/ceed.h>
90d0321e0SJeremy L Thompson #include <ceed/backend.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
110d0321e0SJeremy L Thompson #include <iostream>
120d0321e0SJeremy L Thompson #include <sstream>
130d0321e0SJeremy L Thompson #include <string.h>
140d0321e0SJeremy L Thompson #include "ceed-hip-ref.h"
150d0321e0SJeremy L Thompson #include "../hip/ceed-hip-compile.h"
160d0321e0SJeremy L Thompson 
170d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
180d0321e0SJeremy L Thompson // Build QFunction kernel
190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
200d0321e0SJeremy L Thompson extern "C" int CeedHipBuildQFunction(CeedQFunction qf) {
210d0321e0SJeremy L Thompson   CeedInt ierr;
220d0321e0SJeremy L Thompson   using std::ostringstream;
230d0321e0SJeremy L Thompson   using std::string;
24437930d1SJeremy L Thompson   Ceed ceed;
25437930d1SJeremy L Thompson   CeedQFunctionGetCeed(qf, &ceed);
2686e1ed65Snbeams   Ceed_Hip *ceed_Hip;
2786e1ed65Snbeams   ierr = CeedGetData(ceed, &ceed_Hip); CeedChkBackend(ierr);
280d0321e0SJeremy L Thompson   CeedQFunction_Hip *data;
290d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetData(qf, (void **)&data); CeedChkBackend(ierr);
30437930d1SJeremy L Thompson 
310d0321e0SJeremy L Thompson   // QFunction is built
32437930d1SJeremy L Thompson   if (data->QFunction)
330d0321e0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
340d0321e0SJeremy L Thompson 
35437930d1SJeremy L Thompson   if (!data->qfunction_source)
36437930d1SJeremy L Thompson     // LCOV_EXCL_START
37437930d1SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
38437930d1SJeremy L Thompson                      "No QFunction source or hipFunction_t provided.");
39437930d1SJeremy L Thompson   // LCOV_EXCL_STOP
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;
44437930d1SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, &num_input_fields, &input_fields,
45437930d1SJeremy L Thompson                                 &num_output_fields, &output_fields);
460d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
470d0321e0SJeremy L Thompson 
480d0321e0SJeremy L Thompson   // Build strings for final kernel
49437930d1SJeremy L Thompson   char *read_write_kernel_path, *read_write_kernel_source;
50ee5a26f2SJeremy L Thompson   ierr = CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-ref-qfunction.h",
51437930d1SJeremy L Thompson                                 &read_write_kernel_path); CeedChkBackend(ierr);
5246dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading QFunction Read/Write Kernel Source -----\n");
53437930d1SJeremy L Thompson   ierr = CeedLoadSourceToBuffer(ceed, read_write_kernel_path, &read_write_kernel_source);
54437930d1SJeremy L Thompson   CeedChkBackend(ierr);
5546dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading QFunction Read/Write Kernel Source Complete! -----\n");
56437930d1SJeremy L Thompson   string qfunction_source(data->qfunction_source);
57437930d1SJeremy L Thompson   string qfunction_name(data->qfunction_name);
58437930d1SJeremy L Thompson   string read_write(read_write_kernel_source);
59*204bfdd7SJeremy L Thompson   string kernel_name = "CeedKernelHipRefQFunction_" + qfunction_name;
600d0321e0SJeremy L Thompson   ostringstream code;
610d0321e0SJeremy L Thompson 
620d0321e0SJeremy L Thompson   // Defintions
63437930d1SJeremy L Thompson   code << read_write;
64437930d1SJeremy L Thompson   code << qfunction_source;
6546dc0734SJeremy L Thompson   code << "\n";
6686e1ed65Snbeams   code << "extern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
6786e1ed65Snbeams   code << "__global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Hip fields) {\n";
680d0321e0SJeremy L Thompson 
690d0321e0SJeremy L Thompson   // Inputs
7046dc0734SJeremy L Thompson   code << "  // Input fields\n";
71437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
72437930d1SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(input_fields[i], &size); CeedChkBackend(ierr);
7346dc0734SJeremy L Thompson     code << "  const CeedInt size_input_" << i << " = " << size << ";\n";
7446dc0734SJeremy L Thompson     code << "  CeedScalar input_" << i << "[size_input_" << i << "];\n";
7546dc0734SJeremy L Thompson   }
7646dc0734SJeremy L Thompson   code << "  const CeedScalar* inputs[" << num_input_fields << "];\n";
7746dc0734SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
7846dc0734SJeremy L Thompson     code << "  inputs[" << i << "] = input_" << i << ";\n";
790d0321e0SJeremy L Thompson   }
80437930d1SJeremy L Thompson   code << "\n";
810d0321e0SJeremy L Thompson 
820d0321e0SJeremy L Thompson   // Outputs
8346dc0734SJeremy L Thompson   code << "  // Output fields\n";
84437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
85437930d1SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(output_fields[i], &size); CeedChkBackend(ierr);
8646dc0734SJeremy L Thompson     code << "  const CeedInt size_output_" << i << " = " << size << ";\n";
8746dc0734SJeremy L Thompson     code << "  CeedScalar output_" << i << "[size_output_" << i << "];\n";
880d0321e0SJeremy L Thompson   }
8946dc0734SJeremy L Thompson   code << "  CeedScalar* outputs[" << num_output_fields << "];\n";
90437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
9146dc0734SJeremy L Thompson     code << "  outputs[" << i << "] = output_" << i << ";\n";
920d0321e0SJeremy L Thompson   }
93437930d1SJeremy L Thompson   code << "\n";
940d0321e0SJeremy L Thompson 
950d0321e0SJeremy L Thompson   // Loop over quadrature points
9646dc0734SJeremy L Thompson   code << "  // Loop over quadrature points\n";
970d0321e0SJeremy L Thompson   code << "  for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n";
980d0321e0SJeremy L Thompson 
990d0321e0SJeremy L Thompson   // Load inputs
10046dc0734SJeremy L Thompson   code << "    // -- Load inputs\n";
101437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
10246dc0734SJeremy L Thompson     code << "    readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n";
1030d0321e0SJeremy L Thompson   }
10446dc0734SJeremy L Thompson   code << "\n";
10546dc0734SJeremy L Thompson 
1060d0321e0SJeremy L Thompson   // QFunction
10746dc0734SJeremy L Thompson   code << "    // -- Call QFunction\n";
10846dc0734SJeremy L Thompson   code << "    " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n";
1090d0321e0SJeremy L Thompson 
1100d0321e0SJeremy L Thompson   // Write outputs
11146dc0734SJeremy L Thompson   code << "    // -- Write outputs\n";
112437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
11346dc0734SJeremy L Thompson     code << "    writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n";
1140d0321e0SJeremy L Thompson   }
1150d0321e0SJeremy L Thompson   code << "  }\n";
1160d0321e0SJeremy L Thompson   code << "}\n";
1170d0321e0SJeremy L Thompson 
1180d0321e0SJeremy L Thompson   // View kernel for debugging
11946dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "Generated QFunction Kernels:\n");
1200d0321e0SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
1210d0321e0SJeremy L Thompson 
1220d0321e0SJeremy L Thompson   // Compile kernel
12386e1ed65Snbeams   ierr = CeedCompileHip(ceed, code.str().c_str(), &data->module,
12486e1ed65Snbeams 		        1, "BLOCK_SIZE", ceed_Hip->opt_block_size);
1250d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
126437930d1SJeremy L Thompson   ierr = CeedGetKernelHip(ceed, data->module, kernel_name.c_str(), &data->QFunction);
1270d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
1280d0321e0SJeremy L Thompson 
1290d0321e0SJeremy L Thompson   // Cleanup
130437930d1SJeremy L Thompson   ierr = CeedFree(&data->qfunction_source); CeedChkBackend(ierr);
131437930d1SJeremy L Thompson   ierr = CeedFree(&read_write_kernel_path); CeedChkBackend(ierr);
132437930d1SJeremy L Thompson   ierr = CeedFree(&read_write_kernel_source); CeedChkBackend(ierr);
133437930d1SJeremy L Thompson 
1340d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1350d0321e0SJeremy L Thompson }
1360d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
137