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