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 = CeedGetInstalledJitPath(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 << "\n#define CEED_QFUNCTION(name) inline __device__ int name\n"; 62 code << "#define CEED_QFUNCTION_HELPER inline __device__\n"; 63 code << "#define CeedPragmaSIMD\n"; 64 code << "#define CEED_ERROR_SUCCESS 0\n"; 65 code << "#define CEED_Q_VLA 1\n\n"; 66 code << "typedef struct { const CeedScalar* inputs[16]; CeedScalar* outputs[16]; } Fields_Cuda;\n"; 67 code << read_write; 68 code << qfunction_source; 69 code << "\n"; 70 code << "extern \"C\" __global__ void " << kernel_name << "(void *ctx, CeedInt Q, Fields_Cuda fields) {\n"; 71 72 // Inputs 73 code << " // Input fields\n"; 74 for (CeedInt i = 0; i < num_input_fields; i++) { 75 ierr = CeedQFunctionFieldGetSize(input_fields[i], &size); CeedChkBackend(ierr); 76 code << " const CeedInt size_input_" << i << " = " << size << ";\n"; 77 code << " CeedScalar input_" << i << "[size_input_" << i << "];\n"; 78 } 79 code << " const CeedScalar* inputs[" << num_input_fields << "];\n"; 80 for (CeedInt i = 0; i < num_input_fields; i++) { 81 code << " inputs[" << i << "] = input_" << i << ";\n"; 82 } 83 code << "\n"; 84 85 // Outputs 86 code << " // Output fields\n"; 87 for (CeedInt i = 0; i < num_output_fields; i++) { 88 ierr = CeedQFunctionFieldGetSize(output_fields[i], &size); CeedChkBackend(ierr); 89 code << " const CeedInt size_output_" << i << " = " << size << ";\n"; 90 code << " CeedScalar output_" << i << "[size_output_" << i << "];\n"; 91 } 92 code << " CeedScalar* outputs[" << num_output_fields << "];\n"; 93 for (CeedInt i = 0; i < num_output_fields; i++) { 94 code << " outputs[" << i << "] = output_" << i << ";\n"; 95 } 96 code << "\n"; 97 98 // Loop over quadrature points 99 code << " // Loop over quadrature points\n"; 100 code << " for (CeedInt q = blockIdx.x * blockDim.x + threadIdx.x; q < Q; q += blockDim.x * gridDim.x) {\n"; 101 102 // Load inputs 103 code << " // -- Load inputs\n"; 104 for (CeedInt i = 0; i < num_input_fields; i++) { 105 code << " readQuads<size_input_" << i << ">(q, Q, fields.inputs[" << i << "], input_" << i << ");\n"; 106 } 107 code << "\n"; 108 109 // QFunction 110 code << " // -- Call QFunction\n"; 111 code << " " << qfunction_name << "(ctx, 1, inputs, outputs);\n\n"; 112 113 // Write outputs 114 code << " // -- Write outputs\n"; 115 for (CeedInt i = 0; i < num_output_fields; i++) { 116 code << " writeQuads<size_output_" << i << ">(q, Q, output_" << i << ", fields.outputs[" << i << "]);\n"; 117 } 118 code << " }\n"; 119 code << "}\n"; 120 121 // View kernel for debugging 122 CeedDebug256(ceed, 2, "Generated QFunction Kernels:\n"); 123 CeedDebug(ceed, code.str().c_str()); 124 125 // Compile kernel 126 ierr = CeedCompileCuda(ceed, code.str().c_str(), &data->module, 0); 127 CeedChkBackend(ierr); 128 ierr = CeedGetKernelCuda(ceed, data->module, kernel_name.c_str(), &data->QFunction); 129 CeedChkBackend(ierr); 130 131 // Cleanup 132 ierr = CeedFree(&data->qfunction_source); CeedChkBackend(ierr); 133 ierr = CeedFree(&read_write_kernel_path); CeedChkBackend(ierr); 134 ierr = CeedFree(&read_write_kernel_source); CeedChkBackend(ierr); 135 136 return CEED_ERROR_SUCCESS; 137 } 138 //------------------------------------------------------------------------------ 139