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