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