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