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