xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-qfunction-load.cpp (revision bdee0278611904727ee35fcc2d0d7c3bf83db4c4)
1 // Copyright (c) 2017-2026, 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 CeedQFunctionBuildKernel_Cuda_ref(CeedQFunction qf) {
23   using std::ostringstream;
24   using std::string;
25 
26   Ceed                ceed;
27   CeedInt             num_input_fields, num_output_fields, size;
28   CeedQFunctionField *input_fields, *output_fields;
29   CeedQFunction_Cuda *data;
30 
31   // QFunction is built
32   CeedCallBackend(CeedQFunctionGetData(qf, (void **)&data));
33   if (data->QFunction) return CEED_ERROR_SUCCESS;
34 
35   CeedCallBackend(CeedQFunctionGetCeed(qf, &ceed));
36 
37   // QFunction kernel generation
38   CeedCallBackend(CeedQFunctionGetFields(qf, &num_input_fields, &input_fields, &num_output_fields, &output_fields));
39 
40   // Build strings for final kernel
41   string        qfunction_name(data->qfunction_name);
42   string        kernel_name = "CeedKernelCudaRefQFunction_" + qfunction_name;
43   ostringstream code;
44 
45   // Definitions
46   code << "// QFunction source\n";
47   code << "#include <ceed/jit-source/cuda/cuda-ref-qfunction.h>\n\n";
48   {
49     const char *source_path;
50 
51     CeedCallBackend(CeedQFunctionGetSourcePath(qf, &source_path));
52     CeedCheck(source_path, ceed, CEED_ERROR_BACKEND, "No QFunction source or CUfunction provided.");
53 
54     code << "// User QFunction source\n";
55     code << "#include \"" << source_path << "\"\n\n";
56   }
57   code << "extern \"C\" __global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Cuda fields) {\n";
58 
59   // Inputs
60   code << "  // Input fields\n";
61   for (CeedInt i = 0; i < num_input_fields; i++) {
62     CeedCallBackend(CeedQFunctionFieldGetSize(input_fields[i], &size));
63     code << "  const CeedInt size_input_" << i << " = " << size << ";\n";
64     code << "  CeedScalar input_" << i << "[size_input_" << i << "];\n";
65   }
66   code << "  const CeedScalar *inputs[" << CeedIntMax(num_input_fields, 1) << "];\n";
67   for (CeedInt i = 0; i < num_input_fields; i++) {
68     code << "  inputs[" << i << "] = input_" << i << ";\n";
69   }
70   code << "\n";
71 
72   // Outputs
73   code << "  // Output fields\n";
74   for (CeedInt i = 0; i < num_output_fields; i++) {
75     CeedCallBackend(CeedQFunctionFieldGetSize(output_fields[i], &size));
76     code << "  const CeedInt size_output_" << i << " = " << size << ";\n";
77     code << "  CeedScalar output_" << i << "[size_output_" << i << "];\n";
78   }
79   code << "  CeedScalar *outputs[" << CeedIntMax(num_output_fields, 1) << "];\n";
80   for (CeedInt i = 0; i < num_output_fields; i++) {
81     code << "  outputs[" << i << "] = output_" << i << ";\n";
82   }
83   code << "\n";
84 
85   // Loop over quadrature points
86   code << "  // Loop over quadrature points\n";
87   code << "  for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n";
88 
89   // Load inputs
90   code << "    // -- Load inputs\n";
91   for (CeedInt i = 0; i < num_input_fields; i++) {
92     code << "    readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n";
93   }
94   code << "\n";
95 
96   // QFunction
97   code << "    // -- Call QFunction\n";
98   code << "    " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n";
99 
100   // Write outputs
101   code << "    // -- Write outputs\n";
102   for (CeedInt i = 0; i < num_output_fields; i++) {
103     code << "    writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n";
104   }
105   code << "  }\n";
106   code << "}\n";
107 
108   // Compile kernel
109   CeedCallBackend(CeedCompile_Cuda(ceed, code.str().c_str(), &data->module, 0));
110   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, kernel_name.c_str(), &data->QFunction));
111   CeedCallBackend(CeedDestroy(&ceed));
112   return CEED_ERROR_SUCCESS;
113 }
114 
115 //------------------------------------------------------------------------------
116