xref: /libCEED/backends/hip-ref/ceed-hip-ref-qfunction-load.cpp (revision 257916f8766357eb38205c54507de8baf75f531a)
1 // Copyright (c) 2017-2022, 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/ceed.h>
9 #include <ceed/backend.h>
10 #include <ceed/jit-tools.h>
11 #include <iostream>
12 #include <sstream>
13 #include <string.h>
14 #include "ceed-hip-ref.h"
15 #include "../hip/ceed-hip-compile.h"
16 
17 //------------------------------------------------------------------------------
18 // Build QFunction kernel
19 //------------------------------------------------------------------------------
20 extern "C" int CeedHipBuildQFunction(CeedQFunction qf) {
21   CeedInt ierr;
22   using std::ostringstream;
23   using std::string;
24   Ceed ceed;
25   CeedQFunctionGetCeed(qf, &ceed);
26   Ceed_Hip *ceed_Hip;
27   ierr = CeedGetData(ceed, &ceed_Hip); CeedChkBackend(ierr);
28   CeedQFunction_Hip *data;
29   ierr = CeedQFunctionGetData(qf, (void **)&data); CeedChkBackend(ierr);
30 
31   // QFunction is built
32   if (data->QFunction)
33     return CEED_ERROR_SUCCESS;
34 
35   if (!data->qfunction_source)
36     // LCOV_EXCL_START
37     return CeedError(ceed, CEED_ERROR_BACKEND,
38                      "No QFunction source or hipFunction_t provided.");
39   // LCOV_EXCL_STOP
40 
41   // QFunction kernel generation
42   CeedInt num_input_fields, num_output_fields, size;
43   CeedQFunctionField *input_fields, *output_fields;
44   ierr = CeedQFunctionGetFields(qf, &num_input_fields, &input_fields,
45                                 &num_output_fields, &output_fields);
46   CeedChkBackend(ierr);
47 
48   // Build strings for final kernel
49   char *read_write_kernel_path, *read_write_kernel_source;
50   ierr = CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-ref-qfunction.h",
51                                 &read_write_kernel_path); CeedChkBackend(ierr);
52   CeedDebug256(ceed, 2, "----- Loading QFunction Read/Write Kernel Source -----\n");
53   ierr = CeedLoadSourceToBuffer(ceed, read_write_kernel_path, &read_write_kernel_source);
54   CeedChkBackend(ierr);
55   CeedDebug256(ceed, 2, "----- Loading QFunction Read/Write Kernel Source Complete! -----\n");
56   string qfunction_source(data->qfunction_source);
57   string qfunction_name(data->qfunction_name);
58   string read_write(read_write_kernel_source);
59   string kernel_name = "CeedKernel_Hip_ref_" + qfunction_name;
60   ostringstream code;
61 
62   // Defintions
63   code << read_write;
64   code << qfunction_source;
65   code << "\n";
66   code << "extern \"C\" __launch_bounds__(BLOCK_SIZE)\n";
67   code << "__global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Hip fields) {\n";
68 
69   // Inputs
70   code << "  // Input fields\n";
71   for (CeedInt i = 0; i < num_input_fields; i++) {
72     ierr = CeedQFunctionFieldGetSize(input_fields[i], &size); CeedChkBackend(ierr);
73     code << "  const CeedInt size_input_" << i << " = " << size << ";\n";
74     code << "  CeedScalar input_" << i << "[size_input_" << i << "];\n";
75   }
76   code << "  const CeedScalar* inputs[" << num_input_fields << "];\n";
77   for (CeedInt i = 0; i < num_input_fields; i++) {
78     code << "  inputs[" << i << "] = input_" << i << ";\n";
79   }
80   code << "\n";
81 
82   // Outputs
83   code << "  // Output fields\n";
84   for (CeedInt i = 0; i < num_output_fields; i++) {
85     ierr = CeedQFunctionFieldGetSize(output_fields[i], &size); CeedChkBackend(ierr);
86     code << "  const CeedInt size_output_" << i << " = " << size << ";\n";
87     code << "  CeedScalar output_" << i << "[size_output_" << i << "];\n";
88   }
89   code << "  CeedScalar* outputs[" << num_output_fields << "];\n";
90   for (CeedInt i = 0; i < num_output_fields; i++) {
91     code << "  outputs[" << i << "] = output_" << i << ";\n";
92   }
93   code << "\n";
94 
95   // Loop over quadrature points
96   code << "  // Loop over quadrature points\n";
97   code << "  for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n";
98 
99   // Load inputs
100   code << "    // -- Load inputs\n";
101   for (CeedInt i = 0; i < num_input_fields; i++) {
102     code << "    readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n";
103   }
104   code << "\n";
105 
106   // QFunction
107   code << "    // -- Call QFunction\n";
108   code << "    " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n";
109 
110   // Write outputs
111   code << "    // -- Write outputs\n";
112   for (CeedInt i = 0; i < num_output_fields; i++) {
113     code << "    writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n";
114   }
115   code << "  }\n";
116   code << "}\n";
117 
118   // View kernel for debugging
119   CeedDebug256(ceed, 2, "Generated QFunction Kernels:\n");
120   CeedDebug(ceed, code.str().c_str());
121 
122   // Compile kernel
123   ierr = CeedCompileHip(ceed, code.str().c_str(), &data->module,
124 		        1, "BLOCK_SIZE", ceed_Hip->opt_block_size);
125   CeedChkBackend(ierr);
126   ierr = CeedGetKernelHip(ceed, data->module, kernel_name.c_str(), &data->QFunction);
127   CeedChkBackend(ierr);
128 
129   // Cleanup
130   ierr = CeedFree(&data->qfunction_source); CeedChkBackend(ierr);
131   ierr = CeedFree(&read_write_kernel_path); CeedChkBackend(ierr);
132   ierr = CeedFree(&read_write_kernel_source); CeedChkBackend(ierr);
133 
134   return CEED_ERROR_SUCCESS;
135 }
136 //------------------------------------------------------------------------------
137