xref: /libCEED/rust/libceed-sys/c-src/backends/cuda-gen/ceed-cuda-gen-operator.c (revision b2165e7a2e371018feedcb47974a027ed47e0487)
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>
113d576824SJeremy L Thompson #include <stddef.h>
122b730f8bSJeremy L Thompson 
1349aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h"
146d69246aSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h"
152b730f8bSJeremy L Thompson #include "ceed-cuda-gen-operator-build.h"
162b730f8bSJeremy L Thompson #include "ceed-cuda-gen.h"
17241a4b83SYohann 
18ab213215SJeremy L Thompson //------------------------------------------------------------------------------
19ab213215SJeremy L Thompson // Destroy operator
20ab213215SJeremy L Thompson //------------------------------------------------------------------------------
21241a4b83SYohann static int CeedOperatorDestroy_Cuda_gen(CeedOperator op) {
22241a4b83SYohann   CeedOperator_Cuda_gen *impl;
232b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &impl));
242b730f8bSJeremy L Thompson   CeedCallBackend(CeedFree(&impl));
25e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
26241a4b83SYohann }
27241a4b83SYohann 
282b730f8bSJeremy L Thompson static int Waste(int threads_per_sm, int warp_size, int threads_per_elem, int elems_per_block) {
2939532cebSJed Brown   int useful_threads_per_block = threads_per_elem * elems_per_block;
3039532cebSJed Brown   // round up to nearest multiple of warp_size
31*b2165e7aSSebastian Grimberg   int block_size    = CeedDivUpInt(useful_threads_per_block, warp_size) * warp_size;
3239532cebSJed Brown   int blocks_per_sm = threads_per_sm / block_size;
3339532cebSJed Brown   return threads_per_sm - useful_threads_per_block * blocks_per_sm;
3439532cebSJed Brown }
3539532cebSJed Brown 
36ea61e9acSJeremy L Thompson // Choose the least wasteful block size constrained by blocks_per_sm of max_threads_per_block.
3739532cebSJed Brown //
38ea61e9acSJeremy L Thompson // The x and y part of block[] contains per-element sizes (specified on input) while the z part is number of elements.
3939532cebSJed Brown //
40ea61e9acSJeremy L Thompson // Problem setting: we'd like to make occupancy high with relatively few inactive threads. CUDA (cuOccupancyMaxPotentialBlockSize) can tell us how
4139532cebSJed Brown // many threads can run.
4239532cebSJed Brown //
43ea61e9acSJeremy L Thompson // Note that full occupancy sometimes can't be achieved by one thread block.
44ea61e9acSJeremy L Thompson // For example, an SM might support 1536 threads in total, but only 1024 within a single thread block.
45ea61e9acSJeremy 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
46ea61e9acSJeremy 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
47ea61e9acSJeremy L Thompson // relatively low occupancy and smaller thread blocks, but we solve a reasonably general problem here. Empirically, we find that blocks bigger than
48ea61e9acSJeremy L Thompson // about 256 have higher latency and worse load balancing when the number of elements is modest.
4939532cebSJed Brown //
50ea61e9acSJeremy 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).
51ea61e9acSJeremy L Thompson // They also have a lot of __syncthreads(), which is another point against excessively large thread blocks.
52ea61e9acSJeremy L Thompson // Suppose I have elements with 7x7x7 quadrature points.
53ea61e9acSJeremy L Thompson // This will loop over the last dimension, so we have 7*7=49 threads per element.
54ea61e9acSJeremy L Thompson // Suppose we have two elements = 2*49=98 useful threads.
55ea61e9acSJeremy L Thompson // CUDA schedules in units of full warps (32 threads), so 128 CUDA hardware threads are effectively committed to that block.
56ea61e9acSJeremy L Thompson // Now suppose cuOccupancyMaxPotentialBlockSize returned 352.
57ea61e9acSJeremy 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
58ea61e9acSJeremy L Thompson // hardware threads).
5939532cebSJed Brown //
60ea61e9acSJeremy 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.
61ea61e9acSJeremy L Thompson // Alternatively, we could pack a single block of 7 elements (2*49=343 useful threads) into the 354 slots.
62ea61e9acSJeremy L Thompson // The latter has the least "waste", but __syncthreads() over-synchronizes and it might not pay off relative to smaller blocks.
632b730f8bSJeremy 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],
642b730f8bSJeremy L Thompson                               int *grid) {
6539532cebSJed Brown   const int threads_per_sm   = blocks_per_sm * max_threads_per_block;
6639532cebSJed Brown   const int threads_per_elem = block[0] * block[1];
6739532cebSJed Brown   int       elems_per_block  = 1;
6839532cebSJed Brown   int       waste            = Waste(threads_per_sm, warp_size, threads_per_elem, 1);
692b730f8bSJeremy L Thompson   for (int i = 2; i <= CeedIntMin(max_threads_per_block / threads_per_elem, num_elem); i++) {
7039532cebSJed Brown     int i_waste = Waste(threads_per_sm, warp_size, threads_per_elem, i);
71ea61e9acSJeremy 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
7239532cebSJed Brown     // waste as a smaller one, go ahead and prefer the smaller block.
7339532cebSJed Brown     if (i_waste < waste || (i_waste == waste && threads_per_elem * i <= 128)) {
7439532cebSJed Brown       elems_per_block = i;
7539532cebSJed Brown       waste           = i_waste;
7639532cebSJed Brown     }
7739532cebSJed Brown   }
78ea61e9acSJeremy 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
79ea61e9acSJeremy L Thompson   // check before setting the z-dimension size of the block.
8013516544Snbeams   block[2] = CeedIntMin(elems_per_block, max_threads_z);
81*b2165e7aSSebastian Grimberg   *grid    = CeedDivUpInt(num_elem, elems_per_block);
8239532cebSJed Brown   return CEED_ERROR_SUCCESS;
8339532cebSJed Brown }
8439532cebSJed Brown 
85ea61e9acSJeremy L Thompson // callback for cuOccupancyMaxPotentialBlockSize, providing the amount of dynamic shared memory required for a thread block of size threads.
8639532cebSJed Brown static size_t dynamicSMemSize(int threads) { return threads * sizeof(CeedScalar); }
8739532cebSJed Brown 
88ab213215SJeremy L Thompson //------------------------------------------------------------------------------
89ab213215SJeremy L Thompson // Apply and add to output
90ab213215SJeremy L Thompson //------------------------------------------------------------------------------
912b730f8bSJeremy L Thompson static int CeedOperatorApplyAdd_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) {
92241a4b83SYohann   Ceed ceed;
932b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
9439532cebSJed Brown   Ceed_Cuda *cuda_data;
952b730f8bSJeremy L Thompson   CeedCallBackend(CeedGetData(ceed, &cuda_data));
96241a4b83SYohann   CeedOperator_Cuda_gen *data;
972b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetData(op, &data));
98241a4b83SYohann   CeedQFunction           qf;
99241a4b83SYohann   CeedQFunction_Cuda_gen *qf_data;
1002b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetQFunction(op, &qf));
1012b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetData(qf, &qf_data));
1029e201c85SYohann   CeedInt num_elem, num_input_fields, num_output_fields;
1032b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem));
1049e201c85SYohann   CeedOperatorField *op_input_fields, *op_output_fields;
1052b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields));
1069e201c85SYohann   CeedQFunctionField *qf_input_fields, *qf_output_fields;
1072b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields));
1089e201c85SYohann   CeedEvalMode eval_mode;
1092c2ea1dbSJeremy L Thompson   CeedVector   vec, output_vecs[CEED_FIELD_MAX] = {NULL};
110241a4b83SYohann 
111241a4b83SYohann   // Creation of the operator
112eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op));
113241a4b83SYohann 
114241a4b83SYohann   // Input vectors
1159e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1162b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1179e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1189e201c85SYohann       data->fields.inputs[i] = NULL;
119241a4b83SYohann     } else {
120241a4b83SYohann       // Get input vector
1212b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
1229e201c85SYohann       if (vec == CEED_VECTOR_ACTIVE) vec = input_vec;
1232b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i]));
124241a4b83SYohann     }
125241a4b83SYohann   }
126241a4b83SYohann 
127241a4b83SYohann   // Output vectors
1289e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1292b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1309e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
1319e201c85SYohann       data->fields.outputs[i] = NULL;
132241a4b83SYohann     } else {
133241a4b83SYohann       // Get output vector
1342b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
1359e201c85SYohann       if (vec == CEED_VECTOR_ACTIVE) vec = output_vec;
1369e201c85SYohann       output_vecs[i] = vec;
1373b2939feSjeremylt       // Check for multiple output modes
1383b2939feSjeremylt       CeedInt index = -1;
1393b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
1409e201c85SYohann         if (vec == output_vecs[j]) {
1413b2939feSjeremylt           index = j;
1423b2939feSjeremylt           break;
1433b2939feSjeremylt         }
1443b2939feSjeremylt       }
1453b2939feSjeremylt       if (index == -1) {
1462b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorGetArray(vec, CEED_MEM_DEVICE, &data->fields.outputs[i]));
1473b2939feSjeremylt       } else {
1489e201c85SYohann         data->fields.outputs[i] = data->fields.outputs[index];
1493b2939feSjeremylt       }
150241a4b83SYohann     }
151241a4b83SYohann   }
152241a4b83SYohann 
153777ff853SJeremy L Thompson   // Get context data
1542b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c));
155241a4b83SYohann 
156241a4b83SYohann   // Apply operator
1572b730f8bSJeremy L Thompson   void         *opargs[]  = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W};
158241a4b83SYohann   const CeedInt dim       = data->dim;
1599e201c85SYohann   const CeedInt Q_1d      = data->Q_1d;
1609e201c85SYohann   const CeedInt P_1d      = data->max_P_1d;
1619e201c85SYohann   const CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
16239532cebSJed Brown   int           max_threads_per_block, min_grid_size;
1632b730f8bSJeremy L Thompson   CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000));
1642b730f8bSJeremy L Thompson   int block[3] =
1652b730f8bSJeremy L Thompson       {
1662b730f8bSJeremy L Thompson           thread_1d,
1672b730f8bSJeremy L Thompson           dim < 2 ? 1 : thread_1d,
1682b730f8bSJeremy L Thompson           -1,
1692b730f8bSJeremy L Thompson       },
1702b730f8bSJeremy L Thompson       grid;
1712b730f8bSJeremy L Thompson   CeedChkBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, max_threads_per_block,
1722b730f8bSJeremy L Thompson                                     cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid));
17339532cebSJed Brown   CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar);
174eb7e6cafSJeremy L Thompson   CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->op, grid, block[0], block[1], block[2], shared_mem, opargs));
175241a4b83SYohann 
176241a4b83SYohann   // Restore input arrays
1779e201c85SYohann   for (CeedInt i = 0; i < num_input_fields; i++) {
1782b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode));
1799e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
180241a4b83SYohann     } else {
1812b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec));
1829e201c85SYohann       if (vec == CEED_VECTOR_ACTIVE) vec = input_vec;
1832b730f8bSJeremy L Thompson       CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i]));
184241a4b83SYohann     }
185241a4b83SYohann   }
186241a4b83SYohann 
187241a4b83SYohann   // Restore output arrays
1889e201c85SYohann   for (CeedInt i = 0; i < num_output_fields; i++) {
1892b730f8bSJeremy L Thompson     CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode));
1909e201c85SYohann     if (eval_mode == CEED_EVAL_WEIGHT) {  // Skip
191241a4b83SYohann     } else {
1922b730f8bSJeremy L Thompson       CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec));
1939e201c85SYohann       if (vec == CEED_VECTOR_ACTIVE) vec = output_vec;
1943b2939feSjeremylt       // Check for multiple output modes
1953b2939feSjeremylt       CeedInt index = -1;
1963b2939feSjeremylt       for (CeedInt j = 0; j < i; j++) {
1979e201c85SYohann         if (vec == output_vecs[j]) {
1983b2939feSjeremylt           index = j;
1993b2939feSjeremylt           break;
2003b2939feSjeremylt         }
2013b2939feSjeremylt       }
2023b2939feSjeremylt       if (index == -1) {
2032b730f8bSJeremy L Thompson         CeedCallBackend(CeedVectorRestoreArray(vec, &data->fields.outputs[i]));
204241a4b83SYohann       }
205241a4b83SYohann     }
2063b2939feSjeremylt   }
207777ff853SJeremy L Thompson 
208777ff853SJeremy L Thompson   // Restore context data
2092b730f8bSJeremy L Thompson   CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c));
210441428dfSJeremy L Thompson 
211e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
212241a4b83SYohann }
213241a4b83SYohann 
214ab213215SJeremy L Thompson //------------------------------------------------------------------------------
215ab213215SJeremy L Thompson // Create operator
216ab213215SJeremy L Thompson //------------------------------------------------------------------------------
217241a4b83SYohann int CeedOperatorCreate_Cuda_gen(CeedOperator op) {
218241a4b83SYohann   Ceed ceed;
2192b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorGetCeed(op, &ceed));
220241a4b83SYohann   CeedOperator_Cuda_gen *impl;
221241a4b83SYohann 
2222b730f8bSJeremy L Thompson   CeedCallBackend(CeedCalloc(1, &impl));
2232b730f8bSJeremy L Thompson   CeedCallBackend(CeedOperatorSetData(op, impl));
224241a4b83SYohann 
2252b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda_gen));
2262b730f8bSJeremy L Thompson   CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda_gen));
227e15f9bd0SJeremy L Thompson   return CEED_ERROR_SUCCESS;
228241a4b83SYohann }
2296aa95790SJeremy L Thompson 
230ab213215SJeremy L Thompson //------------------------------------------------------------------------------
231