xref: /libCEED/backends/cuda-ref/ceed-cuda-ref-qfunction-load.cpp (revision ee5a26f213b810c411fd8b16edb331d78b70003c)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
30d0321e0SJeremy L Thompson //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
50d0321e0SJeremy L Thompson //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
70d0321e0SJeremy L Thompson 
80d0321e0SJeremy L Thompson #include <ceed/ceed.h>
90d0321e0SJeremy L Thompson #include <ceed/backend.h>
10437930d1SJeremy L Thompson #include <ceed/jit-tools.h>
110d0321e0SJeremy L Thompson #include <iostream>
120d0321e0SJeremy L Thompson #include <sstream>
130d0321e0SJeremy L Thompson #include <string.h>
140d0321e0SJeremy L Thompson #include "ceed-cuda-ref.h"
150d0321e0SJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
160d0321e0SJeremy L Thompson 
170d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
180d0321e0SJeremy L Thompson // Build QFunction kernel
190d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
200d0321e0SJeremy L Thompson extern "C" int CeedCudaBuildQFunction(CeedQFunction qf) {
210d0321e0SJeremy L Thompson   CeedInt ierr;
220d0321e0SJeremy L Thompson   using std::ostringstream;
230d0321e0SJeremy L Thompson   using std::string;
240d0321e0SJeremy L Thompson   Ceed ceed;
250d0321e0SJeremy L Thompson   CeedQFunctionGetCeed(qf, &ceed);
260d0321e0SJeremy L Thompson   CeedQFunction_Cuda *data;
270d0321e0SJeremy L Thompson   ierr = CeedQFunctionGetData(qf, (void **)&data); CeedChkBackend(ierr);
28437930d1SJeremy L Thompson 
290d0321e0SJeremy L Thompson   // QFunction is built
30437930d1SJeremy L Thompson   if (data->QFunction)
310d0321e0SJeremy L Thompson     return CEED_ERROR_SUCCESS;
32437930d1SJeremy L Thompson 
33437930d1SJeremy L Thompson   if (!data->qfunction_source)
34437930d1SJeremy L Thompson     // LCOV_EXCL_START
35437930d1SJeremy L Thompson     return CeedError(ceed, CEED_ERROR_BACKEND,
36437930d1SJeremy L Thompson                      "No QFunction source or CUfunction provided.");
37437930d1SJeremy L Thompson   // LCOV_EXCL_STOP
380d0321e0SJeremy L Thompson 
390d0321e0SJeremy L Thompson   // QFunction kernel generation
40437930d1SJeremy L Thompson   CeedInt num_input_fields, num_output_fields, size;
41437930d1SJeremy L Thompson   CeedQFunctionField *input_fields, *output_fields;
42437930d1SJeremy L Thompson   ierr = CeedQFunctionGetFields(qf, &num_input_fields, &input_fields,
43437930d1SJeremy L Thompson                                 &num_output_fields, &output_fields);
440d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
450d0321e0SJeremy L Thompson 
460d0321e0SJeremy L Thompson   // Build strings for final kernel
47437930d1SJeremy L Thompson   char *read_write_kernel_path, *read_write_kernel_source;
48*ee5a26f2SJeremy L Thompson   ierr = CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-qfunction.h",
49437930d1SJeremy L Thompson                                 &read_write_kernel_path); CeedChkBackend(ierr);
5046dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading QFunction Read/Write Kernel Source -----\n");
51437930d1SJeremy L Thompson   ierr = CeedLoadSourceToBuffer(ceed, read_write_kernel_path, &read_write_kernel_source);
52437930d1SJeremy L Thompson   CeedChkBackend(ierr);
5346dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "----- Loading QFunction Read/Write Kernel Source Complete! -----\n");
54437930d1SJeremy L Thompson   string qfunction_source(data->qfunction_source);
55437930d1SJeremy L Thompson   string qfunction_name(data->qfunction_name);
56437930d1SJeremy L Thompson   string read_write(read_write_kernel_source);
57437930d1SJeremy L Thompson   string kernel_name = "CeedKernel_Cuda_ref_" + qfunction_name;
580d0321e0SJeremy L Thompson   ostringstream code;
590d0321e0SJeremy L Thompson 
600d0321e0SJeremy L Thompson   // Defintions
610d0321e0SJeremy L Thompson   code << "\n#define CEED_QFUNCTION(name) inline __device__ int name\n";
620d0321e0SJeremy L Thompson   code << "#define CEED_QFUNCTION_HELPER inline __device__\n";
630d0321e0SJeremy L Thompson   code << "#define CeedPragmaSIMD\n";
640d0321e0SJeremy L Thompson   code << "#define CEED_ERROR_SUCCESS 0\n";
650d0321e0SJeremy L Thompson   code << "#define CEED_Q_VLA 1\n\n";
660d0321e0SJeremy L Thompson   code << "typedef struct { const CeedScalar* inputs[16]; CeedScalar* outputs[16]; } Fields_Cuda;\n";
67437930d1SJeremy L Thompson   code << read_write;
68437930d1SJeremy L Thompson   code << qfunction_source;
6946dc0734SJeremy L Thompson   code << "\n";
70437930d1SJeremy L Thompson   code << "extern \"C\" __global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Cuda fields) {\n";
710d0321e0SJeremy L Thompson 
720d0321e0SJeremy L Thompson   // Inputs
7346dc0734SJeremy L Thompson   code << "  // Input fields\n";
74437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
75437930d1SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(input_fields[i], &size); CeedChkBackend(ierr);
7646dc0734SJeremy L Thompson     code << "  const CeedInt size_input_" << i << " = " << size << ";\n";
7746dc0734SJeremy L Thompson     code << "  CeedScalar input_" << i << "[size_input_" << i << "];\n";
7846dc0734SJeremy L Thompson   }
7946dc0734SJeremy L Thompson   code << "  const CeedScalar* inputs[" << num_input_fields << "];\n";
8046dc0734SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
8146dc0734SJeremy L Thompson     code << "  inputs[" << i << "] = input_" << i << ";\n";
820d0321e0SJeremy L Thompson   }
83437930d1SJeremy L Thompson   code << "\n";
840d0321e0SJeremy L Thompson 
850d0321e0SJeremy L Thompson   // Outputs
8646dc0734SJeremy L Thompson   code << "  // Output fields\n";
87437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
88437930d1SJeremy L Thompson     ierr = CeedQFunctionFieldGetSize(output_fields[i], &size); CeedChkBackend(ierr);
8946dc0734SJeremy L Thompson     code << "  const CeedInt size_output_" << i << " = " << size << ";\n";
9046dc0734SJeremy L Thompson     code << "  CeedScalar output_" << i << "[size_output_" << i << "];\n";
910d0321e0SJeremy L Thompson   }
9246dc0734SJeremy L Thompson   code << "  CeedScalar* outputs[" << num_output_fields << "];\n";
93437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
9446dc0734SJeremy L Thompson     code << "  outputs[" << i << "] = output_" << i << ";\n";
950d0321e0SJeremy L Thompson   }
96437930d1SJeremy L Thompson   code << "\n";
970d0321e0SJeremy L Thompson 
980d0321e0SJeremy L Thompson   // Loop over quadrature points
9946dc0734SJeremy L Thompson   code << "  // Loop over quadrature points\n";
1000d0321e0SJeremy L Thompson   code << "  for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n";
1010d0321e0SJeremy L Thompson 
1020d0321e0SJeremy L Thompson   // Load inputs
10346dc0734SJeremy L Thompson   code << "    // -- Load inputs\n";
104437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_input_fields; i++) {
10546dc0734SJeremy L Thompson     code << "    readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n";
1060d0321e0SJeremy L Thompson   }
10746dc0734SJeremy L Thompson   code << "\n";
10846dc0734SJeremy L Thompson 
1090d0321e0SJeremy L Thompson   // QFunction
11046dc0734SJeremy L Thompson   code << "    // -- Call QFunction\n";
11146dc0734SJeremy L Thompson   code << "    " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n";
1120d0321e0SJeremy L Thompson 
1130d0321e0SJeremy L Thompson   // Write outputs
11446dc0734SJeremy L Thompson   code << "    // -- Write outputs\n";
115437930d1SJeremy L Thompson   for (CeedInt i = 0; i < num_output_fields; i++) {
11646dc0734SJeremy L Thompson     code << "    writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n";
1170d0321e0SJeremy L Thompson   }
1180d0321e0SJeremy L Thompson   code << "  }\n";
1190d0321e0SJeremy L Thompson   code << "}\n";
1200d0321e0SJeremy L Thompson 
1210d0321e0SJeremy L Thompson   // View kernel for debugging
12246dc0734SJeremy L Thompson   CeedDebug256(ceed, 2, "Generated QFunction Kernels:\n");
1230d0321e0SJeremy L Thompson   CeedDebug(ceed, code.str().c_str());
1240d0321e0SJeremy L Thompson 
1250d0321e0SJeremy L Thompson   // Compile kernel
1260d0321e0SJeremy L Thompson   ierr = CeedCompileCuda(ceed, code.str().c_str(), &data->module, 0);
1270d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
128437930d1SJeremy L Thompson   ierr = CeedGetKernelCuda(ceed, data->module, kernel_name.c_str(), &data->QFunction);
1290d0321e0SJeremy L Thompson   CeedChkBackend(ierr);
1300d0321e0SJeremy L Thompson 
1310d0321e0SJeremy L Thompson   // Cleanup
132437930d1SJeremy L Thompson   ierr = CeedFree(&data->qfunction_source); CeedChkBackend(ierr);
133437930d1SJeremy L Thompson   ierr = CeedFree(&read_write_kernel_path); CeedChkBackend(ierr);
134437930d1SJeremy L Thompson   ierr = CeedFree(&read_write_kernel_source); CeedChkBackend(ierr);
135437930d1SJeremy L Thompson 
1360d0321e0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
1370d0321e0SJeremy L Thompson }
1380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
139