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