xref: /libCEED/backends/hip-ref/ceed-hip-ref-qfunction-load.cpp (revision 8a31047239c7e4b21930903fe948c2d08c208114)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include <ceed.h>
9 #include <ceed/backend.h>
10 #include <ceed/jit-tools.h>
11 #include <string.h>
12 
13 #include <iostream>
14 #include <sstream>
15 
16 #include "../hip/ceed-hip-common.h"
17 #include "../hip/ceed-hip-compile.h"
18 #include "ceed-hip-ref.h"
19 
20 //------------------------------------------------------------------------------
21 // Build QFunction kernel
22 //------------------------------------------------------------------------------
23 extern "C" int CeedQFunctionBuildKernel_Hip_ref(CeedQFunction qf) {
24   using std::ostringstream;
25   using std::string;
26 
27   Ceed                ceed;
28   char               *read_write_kernel_source;
29   const char         *read_write_kernel_path;
30   Ceed_Hip           *ceed_Hip;
31   CeedInt             num_input_fields, num_output_fields, size;
32   CeedQFunctionField *input_fields, *output_fields;
33   CeedQFunction_Hip  *data;
34 
35   CeedCallBackend(CeedQFunctionGetCeed(qf, &ceed));
36   CeedCallBackend(CeedGetData(ceed, &ceed_Hip));
37   CeedCallBackend(CeedQFunctionGetData(qf, (void **)&data));
38 
39   // QFunction is built
40   if (data->QFunction) return CEED_ERROR_SUCCESS;
41 
42   CeedCheck(data->qfunction_source, ceed, CEED_ERROR_BACKEND, "No QFunction source or hipFunction_t provided.");
43 
44   // QFunction kernel generation
45   CeedCallBackend(CeedQFunctionGetFields(qf, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
46 
47   // Build strings for final kernel
48   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-ref-qfunction.h", &read_write_kernel_path));
49   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading QFunction Read/Write Kernel Source -----\n");
50   CeedCallBackend(CeedLoadSourceToBuffer(ceed, read_write_kernel_path, &read_write_kernel_source));
51   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading QFunction Read/Write Kernel Source Complete! -----\n");
52   string        qfunction_source(data->qfunction_source);
53   string        qfunction_name(data->qfunction_name);
54   string        read_write(read_write_kernel_source);
55   string        kernel_name = "CeedKernelHipRefQFunction_" + qfunction_name;
56   ostringstream code;
57 
58   // Defintions
59   code << read_write;
60   code << qfunction_source;
61   code << "\n";
62   code << "extern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
63   code << "__global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Hip fields) {\n";
64 
65   // Inputs
66   code << "  // Input fields\n";
67   for (CeedInt i = 0; i < num_input_fields; i++) {
68     CeedCallBackend(CeedQFunctionFieldGetSize(input_fields[i], &size));
69     code << "  const CeedInt size_input_" << i << " = " << size << ";\n";
70     code << "  CeedScalar input_" << i << "[size_input_" << i << "];\n";
71   }
72   code << "  const CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
73   for (CeedInt i = 0; i < num_input_fields; i++) {
74     code << "  inputs[" << i << "] = input_" << i << ";\n";
75   }
76   code << "\n";
77 
78   // Outputs
79   code << "  // Output fields\n";
80   for (CeedInt i = 0; i < num_output_fields; i++) {
81     CeedCallBackend(CeedQFunctionFieldGetSize(output_fields[i], &size));
82     code << "  const CeedInt size_output_" << i << " = " << size << ";\n";
83     code << "  CeedScalar output_" << i << "[size_output_" << i << "];\n";
84   }
85   code << "  CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
86   for (CeedInt i = 0; i < num_output_fields; i++) {
87     code << "  outputs[" << i << "] = output_" << i << ";\n";
88   }
89   code << "\n";
90 
91   // Loop over quadrature points
92   code << "  // Loop over quadrature points\n";
93   code << "  for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n";
94 
95   // Load inputs
96   code << "    // -- Load inputs\n";
97   for (CeedInt i = 0; i < num_input_fields; i++) {
98     code << "    readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n";
99   }
100   code << "\n";
101 
102   // QFunction
103   code << "    // -- Call QFunction\n";
104   code << "    " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n";
105 
106   // Write outputs
107   code << "    // -- Write outputs\n";
108   for (CeedInt i = 0; i < num_output_fields; i++) {
109     code << "    writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n";
110   }
111   code << "  }\n";
112   code << "}\n";
113 
114   // View kernel for debugging
115   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated QFunction Kernels:\n");
116   CeedDebug(ceed, code.str().c_str());
117 
118   // Compile kernel
119   CeedCallBackend(CeedCompile_Hip(ceed, code.str().c_str(), &data->module, 1, "BLOCK_SIZE", ceed_Hip->opt_block_size));
120   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, kernel_name.c_str(), &data->QFunction));
121 
122   // Cleanup
123   CeedCallBackend(CeedFree(&data->qfunction_source));
124   CeedCallBackend(CeedFree(&read_write_kernel_path));
125   CeedCallBackend(CeedFree(&read_write_kernel_source));
126   return CEED_ERROR_SUCCESS;
127 }
128 
129 //------------------------------------------------------------------------------
130