1 // Copyright (c) 2017-2024, 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-source/cuda/cuda-types.h> 11 #include <cuda.h> 12 #include <cuda_runtime.h> 13 #include <stddef.h> 14 #include <string.h> 15 16 #include "../cuda/ceed-cuda-common.h" 17 #include "../cuda/ceed-cuda-compile.h" 18 #include "ceed-cuda-gen-operator-build.h" 19 #include "ceed-cuda-gen.h" 20 21 //------------------------------------------------------------------------------ 22 // Destroy operator 23 //------------------------------------------------------------------------------ 24 static int CeedOperatorDestroy_Cuda_gen(CeedOperator op) { 25 Ceed ceed; 26 CeedOperator_Cuda_gen *impl; 27 28 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 29 CeedCallBackend(CeedOperatorGetData(op, &impl)); 30 if (impl->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)impl->points.num_per_elem)); 31 CeedCallBackend(CeedFree(&impl)); 32 CeedCallBackend(CeedDestroy(&ceed)); 33 return CEED_ERROR_SUCCESS; 34 } 35 36 static int Waste(int threads_per_sm, int warp_size, int threads_per_elem, int elems_per_block) { 37 int useful_threads_per_block = threads_per_elem * elems_per_block; 38 // round up to nearest multiple of warp_size 39 int block_size = CeedDivUpInt(useful_threads_per_block, warp_size) * warp_size; 40 int blocks_per_sm = threads_per_sm / block_size; 41 return threads_per_sm - useful_threads_per_block * blocks_per_sm; 42 } 43 44 // Choose the least wasteful block size constrained by blocks_per_sm of max_threads_per_block. 45 // 46 // The x and y part of block[] contains per-element sizes (specified on input) while the z part is number of elements. 47 // 48 // Problem setting: we'd like to make occupancy high with relatively few inactive threads. CUDA (cuOccupancyMaxPotentialBlockSize) can tell us how 49 // many threads can run. 50 // 51 // Note that full occupancy sometimes can't be achieved by one thread block. 52 // For example, an SM might support 1536 threads in total, but only 1024 within a single thread block. 53 // So cuOccupancyMaxPotentialBlockSize may suggest a block size of 768 so that two blocks can run, versus one block of 1024 will prevent a second 54 // block from running. The cuda-gen kernels are pretty heavy with lots of instruction-level parallelism (ILP) so we'll generally be okay with 55 // relatively low occupancy and smaller thread blocks, but we solve a reasonably general problem here. Empirically, we find that blocks bigger than 56 // about 256 have higher latency and worse load balancing when the number of elements is modest. 57 // 58 // cuda-gen can't choose block sizes arbitrarily; they need to be a multiple of the number of quadrature points (or number of basis functions). 59 // They also have a lot of __syncthreads(), which is another point against excessively large thread blocks. 60 // Suppose I have elements with 7x7x7 quadrature points. 61 // This will loop over the last dimension, so we have 7*7=49 threads per element. 62 // Suppose we have two elements = 2*49=98 useful threads. 63 // CUDA schedules in units of full warps (32 threads), so 128 CUDA hardware threads are effectively committed to that block. 64 // Now suppose cuOccupancyMaxPotentialBlockSize returned 352. 65 // We can schedule 2 blocks of size 98 (196 useful threads using 256 hardware threads), but not a third block (which would need a total of 384 66 // hardware threads). 67 // 68 // If instead, we had packed 3 elements, we'd have 3*49=147 useful threads occupying 160 slots, and could schedule two blocks. 69 // Alternatively, we could pack a single block of 7 elements (2*49=343 useful threads) into the 354 slots. 70 // The latter has the least "waste", but __syncthreads() over-synchronizes and it might not pay off relative to smaller blocks. 71 static int BlockGridCalculate(CeedInt num_elem, int blocks_per_sm, int max_threads_per_block, int max_threads_z, int warp_size, int block[3], 72 int *grid) { 73 const int threads_per_sm = blocks_per_sm * max_threads_per_block; 74 const int threads_per_elem = block[0] * block[1]; 75 int elems_per_block = 1; 76 int waste = Waste(threads_per_sm, warp_size, threads_per_elem, 1); 77 78 for (int i = 2; i <= CeedIntMin(max_threads_per_block / threads_per_elem, num_elem); i++) { 79 int i_waste = Waste(threads_per_sm, warp_size, threads_per_elem, i); 80 81 // We want to minimize waste, but smaller kernels have lower latency and less __syncthreads() overhead so when a larger block size has the same 82 // waste as a smaller one, go ahead and prefer the smaller block. 83 if (i_waste < waste || (i_waste == waste && threads_per_elem * i <= 128)) { 84 elems_per_block = i; 85 waste = i_waste; 86 } 87 } 88 // In low-order elements, threads_per_elem may be sufficiently low to give an elems_per_block greater than allowable for the device, so we must 89 // check before setting the z-dimension size of the block. 90 block[2] = CeedIntMin(elems_per_block, max_threads_z); 91 *grid = CeedDivUpInt(num_elem, elems_per_block); 92 return CEED_ERROR_SUCCESS; 93 } 94 95 // callback for cuOccupancyMaxPotentialBlockSize, providing the amount of dynamic shared memory required for a thread block of size threads. 96 static size_t dynamicSMemSize(int threads) { return threads * sizeof(CeedScalar); } 97 98 //------------------------------------------------------------------------------ 99 // Apply and add to output 100 //------------------------------------------------------------------------------ 101 static int CeedOperatorApplyAdd_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) { 102 bool is_at_points, is_tensor, is_run_good = true; 103 Ceed ceed; 104 Ceed_Cuda *cuda_data; 105 CeedInt num_elem, num_input_fields, num_output_fields; 106 CeedEvalMode eval_mode; 107 CeedVector output_vecs[CEED_FIELD_MAX] = {NULL}; 108 CeedQFunctionField *qf_input_fields, *qf_output_fields; 109 CeedQFunction_Cuda_gen *qf_data; 110 CeedQFunction qf; 111 CeedOperatorField *op_input_fields, *op_output_fields; 112 CeedOperator_Cuda_gen *data; 113 114 // Creation of the operator 115 { 116 bool is_good_build = false; 117 118 CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op, &is_good_build)); 119 if (!is_good_build) { 120 CeedOperator op_fallback; 121 122 CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator due to code generation issue"); 123 CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback)); 124 CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request)); 125 return CEED_ERROR_SUCCESS; 126 } 127 } 128 129 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 130 CeedCallBackend(CeedGetData(ceed, &cuda_data)); 131 CeedCallBackend(CeedOperatorGetData(op, &data)); 132 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 133 CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 134 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 135 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 136 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 137 138 // Input vectors 139 for (CeedInt i = 0; i < num_input_fields; i++) { 140 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 141 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 142 data->fields.inputs[i] = NULL; 143 } else { 144 bool is_active; 145 CeedVector vec; 146 147 // Get input vector 148 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 149 is_active = vec == CEED_VECTOR_ACTIVE; 150 if (is_active) vec = input_vec; 151 CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i])); 152 if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec)); 153 } 154 } 155 156 // Output vectors 157 for (CeedInt i = 0; i < num_output_fields; i++) { 158 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 159 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 160 data->fields.outputs[i] = NULL; 161 } else { 162 bool is_active; 163 CeedVector vec; 164 165 // Get output vector 166 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 167 is_active = vec == CEED_VECTOR_ACTIVE; 168 if (is_active) vec = output_vec; 169 output_vecs[i] = vec; 170 // Check for multiple output modes 171 CeedInt index = -1; 172 173 for (CeedInt j = 0; j < i; j++) { 174 if (vec == output_vecs[j]) { 175 index = j; 176 break; 177 } 178 } 179 if (index == -1) { 180 CeedCallBackend(CeedVectorGetArray(vec, CEED_MEM_DEVICE, &data->fields.outputs[i])); 181 } else { 182 data->fields.outputs[i] = data->fields.outputs[index]; 183 } 184 if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec)); 185 } 186 } 187 188 // Point coordinates, if needed 189 CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points)); 190 if (is_at_points) { 191 // Coords 192 CeedVector vec; 193 194 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec)); 195 CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords)); 196 CeedCallBackend(CeedVectorDestroy(&vec)); 197 198 // Points per elem 199 if (num_elem != data->points.num_elem) { 200 CeedInt *points_per_elem; 201 const CeedInt num_bytes = num_elem * sizeof(CeedInt); 202 CeedElemRestriction rstr_points = NULL; 203 204 data->points.num_elem = num_elem; 205 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL)); 206 CeedCallBackend(CeedCalloc(num_elem, &points_per_elem)); 207 for (CeedInt e = 0; e < num_elem; e++) { 208 CeedInt num_points_elem; 209 210 CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem)); 211 points_per_elem[e] = num_points_elem; 212 } 213 if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem)); 214 CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes)); 215 CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice)); 216 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points)); 217 CeedCallBackend(CeedFree(&points_per_elem)); 218 } 219 } 220 221 // Get context data 222 CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c)); 223 224 // Apply operator 225 void *opargs[] = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W, &data->points}; 226 const CeedInt dim = data->dim; 227 const CeedInt Q_1d = data->Q_1d; 228 const CeedInt P_1d = data->max_P_1d; 229 const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 230 int max_threads_per_block, min_grid_size, grid; 231 232 CeedCallBackend(CeedOperatorHasTensorBases(op, &is_tensor)); 233 CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000)); 234 int block[3] = {thread_1d, ((!is_tensor || dim == 1) ? 1 : thread_1d), -1}; 235 236 if (is_tensor) { 237 CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, max_threads_per_block, 238 cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid)); 239 } else { 240 CeedInt elems_per_block = CeedIntMin(cuda_data->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 1)); 241 242 grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 243 block[2] = elems_per_block; 244 } 245 CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar); 246 247 CeedCallBackend(CeedTryRunKernelDimShared_Cuda(ceed, data->op, grid, block[0], block[1], block[2], shared_mem, &is_run_good, opargs)); 248 249 // Restore input arrays 250 for (CeedInt i = 0; i < num_input_fields; i++) { 251 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 252 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 253 } else { 254 bool is_active; 255 CeedVector vec; 256 257 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 258 is_active = vec == CEED_VECTOR_ACTIVE; 259 if (is_active) vec = input_vec; 260 CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i])); 261 if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec)); 262 } 263 } 264 265 // Restore output arrays 266 for (CeedInt i = 0; i < num_output_fields; i++) { 267 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 268 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 269 } else { 270 bool is_active; 271 CeedVector vec; 272 273 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 274 is_active = vec == CEED_VECTOR_ACTIVE; 275 if (is_active) vec = output_vec; 276 // Check for multiple output modes 277 CeedInt index = -1; 278 279 for (CeedInt j = 0; j < i; j++) { 280 if (vec == output_vecs[j]) { 281 index = j; 282 break; 283 } 284 } 285 if (index == -1) { 286 CeedCallBackend(CeedVectorRestoreArray(vec, &data->fields.outputs[i])); 287 } 288 if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec)); 289 } 290 } 291 292 // Restore point coordinates, if needed 293 if (is_at_points) { 294 CeedVector vec; 295 296 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec)); 297 CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords)); 298 CeedCallBackend(CeedVectorDestroy(&vec)); 299 } 300 301 // Restore context data 302 CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c)); 303 304 // Cleanup 305 CeedCallBackend(CeedDestroy(&ceed)); 306 CeedCallBackend(CeedQFunctionDestroy(&qf)); 307 308 // Fallback if run was bad (out of resources) 309 if (!is_run_good) { 310 CeedOperator op_fallback; 311 312 data->use_fallback = true; 313 CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator due to kernel execution issue"); 314 CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback)); 315 CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request)); 316 return CEED_ERROR_SUCCESS; 317 } 318 return CEED_ERROR_SUCCESS; 319 } 320 321 //------------------------------------------------------------------------------ 322 // Create operator 323 //------------------------------------------------------------------------------ 324 int CeedOperatorCreate_Cuda_gen(CeedOperator op) { 325 Ceed ceed; 326 CeedOperator_Cuda_gen *impl; 327 328 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 329 CeedCallBackend(CeedCalloc(1, &impl)); 330 CeedCallBackend(CeedOperatorSetData(op, impl)); 331 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda_gen)); 332 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda_gen)); 333 CeedCallBackend(CeedDestroy(&ceed)); 334 return CEED_ERROR_SUCCESS; 335 } 336 337 //------------------------------------------------------------------------------ 338