xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator.c (revision eb7e6cafeb5d6cc94d59355f95e7bc9ae3fc1c25)
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.
3241a4b83SYohann //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5241a4b83SYohann //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7241a4b83SYohann 
849aac155SJeremy L Thompson #include <ceed.h>
9ec3da8bcSJed Brown #include <ceed/backend.h>
1049aac155SJeremy L Thompson #include <ceed/jit-source/cuda/cuda-types.h>
1149aac155SJeremy L Thompson #include <cuda.h>
123d576824SJeremy L Thompson #include <stddef.h>
132b730f8bSJeremy L Thompson 
1449aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
156d69246aSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
162b730f8bSJeremy L Thompson #include "ceed-cuda-gen-operator-build.h"
172b730f8bSJeremy L Thompson #include "ceed-cuda-gen.h"
18241a4b83SYohann 
19ab213215SJeremy L Thompson //------------------------------------------------------------------------------
20ab213215SJeremy L Thompson // Destroy operator
21ab213215SJeremy L Thompson //------------------------------------------------------------------------------
22241a4b83SYohann static int CeedOperatorDestroy_Cuda_gen(CeedOperator op) {
23241a4b83SYohann   CeedOperator_Cuda_gen *impl;
242b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
252b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
26e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
27241a4b83SYohann }
28241a4b83SYohann 
292b730f8bSJeremy L Thompson static int Waste(int threads_per_sm, int warp_size, int threads_per_elem, int elems_per_block) {
3039532cebSJed Brown   int useful_threads_per_block = threads_per_elem * elems_per_block;
3139532cebSJed Brown   // round up to nearest multiple of warp_size
322b730f8bSJeremy L Thompson   int block_size    = ((useful_threads_per_block + warp_size - 1) / warp_size) * warp_size;
3339532cebSJed Brown   int blocks_per_sm = threads_per_sm / block_size;
3439532cebSJed Brown   return threads_per_sm - useful_threads_per_block * blocks_per_sm;
3539532cebSJed Brown }
3639532cebSJed Brown 
37ea61e9acSJeremy L Thompson // Choose the least wasteful block size constrained by blocks_per_sm of max_threads_per_block.
3839532cebSJed Brown //
39ea61e9acSJeremy L Thompson // The x and y part of block[] contains per-element sizes (specified on input) while the z part is number of elements.
4039532cebSJed Brown //
41ea61e9acSJeremy L Thompson // Problem setting: we'd like to make occupancy high with relatively few inactive threads. CUDA (cuOccupancyMaxPotentialBlockSize) can tell us how
4239532cebSJed Brown // many threads can run.
4339532cebSJed Brown //
44ea61e9acSJeremy L Thompson // Note that full occupancy sometimes can't be achieved by one thread block.
45ea61e9acSJeremy L Thompson // For example, an SM might support 1536 threads in total, but only 1024 within a single thread block.
46ea61e9acSJeremy L Thompson // So cuOccupancyMaxPotentialBlockSize may suggest a block size of 768 so that two blocks can run, versus one block of 1024 will prevent a second
47ea61e9acSJeremy L Thompson // block from running. The cuda-gen kernels are pretty heavy with lots of instruction-level parallelism (ILP) so we'll generally be okay with
48ea61e9acSJeremy L Thompson // relatively low occupancy and smaller thread blocks, but we solve a reasonably general problem here. Empirically, we find that blocks bigger than
49ea61e9acSJeremy L Thompson // about 256 have higher latency and worse load balancing when the number of elements is modest.
5039532cebSJed Brown //
51ea61e9acSJeremy L Thompson // 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).
52ea61e9acSJeremy L Thompson // They also have a lot of __syncthreads(), which is another point against excessively large thread blocks.
53ea61e9acSJeremy L Thompson // Suppose I have elements with 7x7x7 quadrature points.
54ea61e9acSJeremy L Thompson // This will loop over the last dimension, so we have 7*7=49 threads per element.
55ea61e9acSJeremy L Thompson // Suppose we have two elements = 2*49=98 useful threads.
56ea61e9acSJeremy L Thompson // CUDA schedules in units of full warps (32 threads), so 128 CUDA hardware threads are effectively committed to that block.
57ea61e9acSJeremy L Thompson // Now suppose cuOccupancyMaxPotentialBlockSize returned 352.
58ea61e9acSJeremy L Thompson // 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
59ea61e9acSJeremy L Thompson // hardware threads).
6039532cebSJed Brown //
61ea61e9acSJeremy L Thompson // If instead, we had packed 3 elements, we'd have 3*49=147 useful threads occupying 160 slots, and could schedule two blocks.
62ea61e9acSJeremy L Thompson // Alternatively, we could pack a single block of 7 elements (2*49=343 useful threads) into the 354 slots.
63ea61e9acSJeremy L Thompson // The latter has the least "waste", but __syncthreads() over-synchronizes and it might not pay off relative to smaller blocks.
642b730f8bSJeremy L Thompson 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],
652b730f8bSJeremy L Thompson                               int *grid) {
6639532cebSJed Brown   const int threads_per_sm   = blocks_per_sm * max_threads_per_block;
6739532cebSJed Brown   const int threads_per_elem = block[0] * block[1];
6839532cebSJed Brown   int       elems_per_block  = 1;
6939532cebSJed Brown   int       waste            = Waste(threads_per_sm, warp_size, threads_per_elem, 1);
702b730f8bSJeremy L Thompson   for (int i = 2; i <= CeedIntMin(max_threads_per_block / threads_per_elem, num_elem); i++) {
7139532cebSJed Brown     int i_waste = Waste(threads_per_sm, warp_size, threads_per_elem, i);
72ea61e9acSJeremy L Thompson     // We want to minimize waste, but smaller kernels have lower latency and less __syncthreads() overhead so when a larger block size has the same
7339532cebSJed Brown     // waste as a smaller one, go ahead and prefer the smaller block.
7439532cebSJed Brown     if (i_waste < waste || (i_waste == waste && threads_per_elem * i <= 128)) {
7539532cebSJed Brown       elems_per_block = i;
7639532cebSJed Brown       waste           = i_waste;
7739532cebSJed Brown     }
7839532cebSJed Brown   }
79ea61e9acSJeremy L Thompson   // 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
80ea61e9acSJeremy L Thompson   // check before setting the z-dimension size of the block.
8113516544Snbeams   block[2] = CeedIntMin(elems_per_block, max_threads_z);
829e201c85SYohann   *grid    = (num_elem + elems_per_block - 1) / elems_per_block;
8339532cebSJed Brown   return CEED_ERROR_SUCCESS;
8439532cebSJed Brown }
8539532cebSJed Brown 
86ea61e9acSJeremy L Thompson // callback for cuOccupancyMaxPotentialBlockSize, providing the amount of dynamic shared memory required for a thread block of size threads.
8739532cebSJed Brown static size_t dynamicSMemSize(int threads) { return threads * sizeof(CeedScalar); }
8839532cebSJed Brown 
89ab213215SJeremy L Thompson //------------------------------------------------------------------------------
90ab213215SJeremy L Thompson // Apply and add to output
91ab213215SJeremy L Thompson //------------------------------------------------------------------------------
922b730f8bSJeremy L Thompson static int CeedOperatorApplyAdd_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) {
93241a4b83SYohann   Ceed ceed;
942b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
9539532cebSJed Brown   Ceed_Cuda *cuda_data;
962b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
97241a4b83SYohann   CeedOperator_Cuda_gen *data;
982b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
99241a4b83SYohann   CeedQFunction           qf;
100241a4b83SYohann   CeedQFunction_Cuda_gen *qf_data;
1012b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1022b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1039e201c85SYohann   CeedInt num_elem, num_input_fields, num_output_fields;
1042b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
1059e201c85SYohann   CeedOperatorField *op_input_fields, *op_output_fields;
1062b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
1079e201c85SYohann   CeedQFunctionField *qf_input_fields, *qf_output_fields;
1082b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1099e201c85SYohann   CeedEvalMode eval_mode;
1102c2ea1dbSJeremy L Thompson   CeedVector   vec, output_vecs[CEED_FIELD_MAX] = {NULL};
111241a4b83SYohann 
112241a4b83SYohann   // Creation of the operator
113*eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op));
114241a4b83SYohann 
115241a4b83SYohann   // Input vectors
1169e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1172b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1189e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1199e201c85SYohann       data->fields.inputs[i] = NULL;
120241a4b83SYohann     } else {
121241a4b83SYohann       // Get input vector
1222b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
1239e201c85SYohann       if (vec == CEED_VECTOR_ACTIVE) vec = input_vec;
1242b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i]));
125241a4b83SYohann     }
126241a4b83SYohann   }
127241a4b83SYohann 
128241a4b83SYohann   // Output vectors
1299e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1302b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1319e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1329e201c85SYohann       data->fields.outputs[i] = NULL;
133241a4b83SYohann     } else {
134241a4b83SYohann       // Get output vector
1352b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
1369e201c85SYohann       if (vec == CEED_VECTOR_ACTIVE) vec = output_vec;
1379e201c85SYohann       output_vecs[i] = vec;
1383b2939feSjeremylt       // Check for multiple output modes
1393b2939feSjeremylt       CeedInt index = -1;
1403b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
1419e201c85SYohann         if (vec == output_vecs[j]) {
1423b2939feSjeremylt           index = j;
1433b2939feSjeremylt           break;
1443b2939feSjeremylt         }
1453b2939feSjeremylt       }
1463b2939feSjeremylt       if (index == -1) {
1472b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(vec, CEED_MEM_DEVICE, &data->fields.outputs[i]));
1483b2939feSjeremylt       } else {
1499e201c85SYohann         data->fields.outputs[i] = data->fields.outputs[index];
1503b2939feSjeremylt       }
151241a4b83SYohann     }
152241a4b83SYohann   }
153241a4b83SYohann 
154777ff853SJeremy L Thompson   // Get context data
1552b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c));
156241a4b83SYohann 
157241a4b83SYohann   // Apply operator
1582b730f8bSJeremy L Thompson 
1592b730f8bSJeremy L Thompson   void         *opargs[]  = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W};
160241a4b83SYohann   const CeedInt dim       = data->dim;
1619e201c85SYohann   const CeedInt Q_1d      = data->Q_1d;
1629e201c85SYohann   const CeedInt P_1d      = data->max_P_1d;
1639e201c85SYohann   const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
16439532cebSJed Brown   int           max_threads_per_block, min_grid_size;
1652b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000));
1662b730f8bSJeremy L Thompson   int block[3] =
1672b730f8bSJeremy L Thompson       {
1682b730f8bSJeremy L Thompson           thread_1d,
1692b730f8bSJeremy L Thompson           dim < 2 ? 1 : thread_1d,
1702b730f8bSJeremy L Thompson           -1,
1712b730f8bSJeremy L Thompson       },
1722b730f8bSJeremy L Thompson       grid;
1732b730f8bSJeremy L Thompson   CeedChkBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, max_threads_per_block,
1742b730f8bSJeremy L Thompson                                     cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid));
17539532cebSJed Brown   CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar);
176*eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->op, grid, block[0], block[1], block[2], shared_mem, opargs));
177241a4b83SYohann 
178241a4b83SYohann   // Restore input arrays
1799e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1802b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1819e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
182241a4b83SYohann     } else {
1832b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
1849e201c85SYohann       if (vec == CEED_VECTOR_ACTIVE) vec = input_vec;
1852b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i]));
186241a4b83SYohann     }
187241a4b83SYohann   }
188241a4b83SYohann 
189241a4b83SYohann   // Restore output arrays
1909e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1912b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1929e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
193241a4b83SYohann     } else {
1942b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
1959e201c85SYohann       if (vec == CEED_VECTOR_ACTIVE) vec = output_vec;
1963b2939feSjeremylt       // Check for multiple output modes
1973b2939feSjeremylt       CeedInt index = -1;
1983b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
1999e201c85SYohann         if (vec == output_vecs[j]) {
2003b2939feSjeremylt           index = j;
2013b2939feSjeremylt           break;
2023b2939feSjeremylt         }
2033b2939feSjeremylt       }
2043b2939feSjeremylt       if (index == -1) {
2052b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(vec, &data->fields.outputs[i]));
206241a4b83SYohann       }
207241a4b83SYohann     }
2083b2939feSjeremylt   }
209777ff853SJeremy L Thompson 
210777ff853SJeremy L Thompson   // Restore context data
2112b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c));
212441428dfSJeremy L Thompson 
213e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
214241a4b83SYohann }
215241a4b83SYohann 
216ab213215SJeremy L Thompson //------------------------------------------------------------------------------
217ab213215SJeremy L Thompson // Create operator
218ab213215SJeremy L Thompson //------------------------------------------------------------------------------
219241a4b83SYohann int CeedOperatorCreate_Cuda_gen(CeedOperator op) {
220241a4b83SYohann   Ceed ceed;
2212b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
222241a4b83SYohann   CeedOperator_Cuda_gen *impl;
223241a4b83SYohann 
2242b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
2252b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
226241a4b83SYohann 
2272b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda_gen));
2282b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda_gen));
229e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
230241a4b83SYohann }
2316aa95790SJeremy L Thompson 
232ab213215SJeremy L Thompson //------------------------------------------------------------------------------
233