1d275d636SJeremy L Thompson // Copyright (c) 2017-2025, 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> 118b97b69aSJeremy L Thompson #include <cuda.h> 128b97b69aSJeremy L Thompson #include <cuda_runtime.h> 133d576824SJeremy L Thompson #include <stddef.h> 14dc007f05SJeremy L Thompson #include <string.h> 152b730f8bSJeremy L Thompson 1649aac155SJeremy L Thompson #include "../cuda/ceed-cuda-common.h" 176d69246aSJeremy L Thompson #include "../cuda/ceed-cuda-compile.h" 182b730f8bSJeremy L Thompson #include "ceed-cuda-gen-operator-build.h" 192b730f8bSJeremy L Thompson #include "ceed-cuda-gen.h" 20241a4b83SYohann 21ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 22ab213215SJeremy L Thompson // Destroy operator 23ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 24241a4b83SYohann static int CeedOperatorDestroy_Cuda_gen(CeedOperator op) { 258b97b69aSJeremy L Thompson Ceed ceed; 26241a4b83SYohann CeedOperator_Cuda_gen *impl; 27ca735530SJeremy L Thompson 288b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 292b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorGetData(op, &impl)); 308b7d3340SJeremy L Thompson if (impl->module) CeedCallCuda(ceed, cuModuleUnload(impl->module)); 31*0183ed61SJeremy L Thompson if (impl->module_assemble_full) CeedCallCuda(ceed, cuModuleUnload(impl->module_assemble_full)); 32*0183ed61SJeremy L Thompson if (impl->module_assemble_diagonal) CeedCallCuda(ceed, cuModuleUnload(impl->module_assemble_diagonal)); 338b97b69aSJeremy L Thompson if (impl->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)impl->points.num_per_elem)); 342b730f8bSJeremy L Thompson CeedCallBackend(CeedFree(&impl)); 358b97b69aSJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 36e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 37241a4b83SYohann } 38241a4b83SYohann 392b730f8bSJeremy L Thompson static int Waste(int threads_per_sm, int warp_size, int threads_per_elem, int elems_per_block) { 4039532cebSJed Brown int useful_threads_per_block = threads_per_elem * elems_per_block; 4139532cebSJed Brown // round up to nearest multiple of warp_size 42b2165e7aSSebastian Grimberg int block_size = CeedDivUpInt(useful_threads_per_block, warp_size) * warp_size; 4339532cebSJed Brown int blocks_per_sm = threads_per_sm / block_size; 4439532cebSJed Brown return threads_per_sm - useful_threads_per_block * blocks_per_sm; 4539532cebSJed Brown } 4639532cebSJed Brown 47ea61e9acSJeremy L Thompson // Choose the least wasteful block size constrained by blocks_per_sm of max_threads_per_block. 4839532cebSJed Brown // 49ea61e9acSJeremy L Thompson // The x and y part of block[] contains per-element sizes (specified on input) while the z part is number of elements. 5039532cebSJed Brown // 51ea61e9acSJeremy L Thompson // Problem setting: we'd like to make occupancy high with relatively few inactive threads. CUDA (cuOccupancyMaxPotentialBlockSize) can tell us how 5239532cebSJed Brown // many threads can run. 5339532cebSJed Brown // 54ea61e9acSJeremy L Thompson // Note that full occupancy sometimes can't be achieved by one thread block. 55ea61e9acSJeremy L Thompson // For example, an SM might support 1536 threads in total, but only 1024 within a single thread block. 56ea61e9acSJeremy 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 57ea61e9acSJeremy 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 58ea61e9acSJeremy L Thompson // relatively low occupancy and smaller thread blocks, but we solve a reasonably general problem here. Empirically, we find that blocks bigger than 59ea61e9acSJeremy L Thompson // about 256 have higher latency and worse load balancing when the number of elements is modest. 6039532cebSJed Brown // 61ea61e9acSJeremy 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). 62ea61e9acSJeremy L Thompson // They also have a lot of __syncthreads(), which is another point against excessively large thread blocks. 63ea61e9acSJeremy L Thompson // Suppose I have elements with 7x7x7 quadrature points. 64ea61e9acSJeremy L Thompson // This will loop over the last dimension, so we have 7*7=49 threads per element. 65ea61e9acSJeremy L Thompson // Suppose we have two elements = 2*49=98 useful threads. 66ea61e9acSJeremy L Thompson // CUDA schedules in units of full warps (32 threads), so 128 CUDA hardware threads are effectively committed to that block. 67ea61e9acSJeremy L Thompson // Now suppose cuOccupancyMaxPotentialBlockSize returned 352. 68ea61e9acSJeremy 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 69ea61e9acSJeremy L Thompson // hardware threads). 7039532cebSJed Brown // 71ea61e9acSJeremy 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. 72ea61e9acSJeremy L Thompson // Alternatively, we could pack a single block of 7 elements (2*49=343 useful threads) into the 354 slots. 73ea61e9acSJeremy L Thompson // The latter has the least "waste", but __syncthreads() over-synchronizes and it might not pay off relative to smaller blocks. 742b730f8bSJeremy 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], 752b730f8bSJeremy L Thompson int *grid) { 7639532cebSJed Brown const int threads_per_sm = blocks_per_sm * max_threads_per_block; 7739532cebSJed Brown const int threads_per_elem = block[0] * block[1]; 7839532cebSJed Brown int elems_per_block = 1; 7939532cebSJed Brown int waste = Waste(threads_per_sm, warp_size, threads_per_elem, 1); 80ca735530SJeremy L Thompson 812b730f8bSJeremy L Thompson for (int i = 2; i <= CeedIntMin(max_threads_per_block / threads_per_elem, num_elem); i++) { 8239532cebSJed Brown int i_waste = Waste(threads_per_sm, warp_size, threads_per_elem, i); 83ca735530SJeremy L Thompson 84ea61e9acSJeremy 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 8539532cebSJed Brown // waste as a smaller one, go ahead and prefer the smaller block. 8639532cebSJed Brown if (i_waste < waste || (i_waste == waste && threads_per_elem * i <= 128)) { 8739532cebSJed Brown elems_per_block = i; 8839532cebSJed Brown waste = i_waste; 8939532cebSJed Brown } 9039532cebSJed Brown } 91ea61e9acSJeremy 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 92ea61e9acSJeremy L Thompson // check before setting the z-dimension size of the block. 9313516544Snbeams block[2] = CeedIntMin(elems_per_block, max_threads_z); 94b2165e7aSSebastian Grimberg *grid = CeedDivUpInt(num_elem, elems_per_block); 9539532cebSJed Brown return CEED_ERROR_SUCCESS; 9639532cebSJed Brown } 9739532cebSJed Brown 98ea61e9acSJeremy L Thompson // callback for cuOccupancyMaxPotentialBlockSize, providing the amount of dynamic shared memory required for a thread block of size threads. 9939532cebSJed Brown static size_t dynamicSMemSize(int threads) { return threads * sizeof(CeedScalar); } 10039532cebSJed Brown 101ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 102ab213215SJeremy L Thompson // Apply and add to output 103ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 104e9c76bddSJeremy L Thompson static int CeedOperatorApplyAddCore_Cuda_gen(CeedOperator op, CUstream stream, const CeedScalar *input_arr, CeedScalar *output_arr, bool *is_run_good, 105ea04d07fSJeremy L Thompson CeedRequest *request) { 106ea04d07fSJeremy L Thompson bool is_at_points, is_tensor; 107241a4b83SYohann Ceed ceed; 10839532cebSJed Brown Ceed_Cuda *cuda_data; 109ca735530SJeremy L Thompson CeedInt num_elem, num_input_fields, num_output_fields; 110ca735530SJeremy L Thompson CeedEvalMode eval_mode; 111ca735530SJeremy L Thompson CeedQFunctionField *qf_input_fields, *qf_output_fields; 112241a4b83SYohann CeedQFunction_Cuda_gen *qf_data; 113ca735530SJeremy L Thompson CeedQFunction qf; 114ca735530SJeremy L Thompson CeedOperatorField *op_input_fields, *op_output_fields; 115ca735530SJeremy L Thompson CeedOperator_Cuda_gen *data; 116ca735530SJeremy L Thompson 117ea04d07fSJeremy L Thompson // Build the operator kernel 118ea04d07fSJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op, is_run_good)); 119ea04d07fSJeremy L Thompson if (!(*is_run_good)) return CEED_ERROR_SUCCESS; 120f6eafd79SJeremy L Thompson 121c11e12f4SJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 122c11e12f4SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &cuda_data)); 123c11e12f4SJeremy L Thompson CeedCallBackend(CeedOperatorGetData(op, &data)); 124c11e12f4SJeremy L Thompson CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 125c11e12f4SJeremy L Thompson CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 126c11e12f4SJeremy L Thompson CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 127ddae5012SJeremy L Thompson CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 128c11e12f4SJeremy L Thompson CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 129c11e12f4SJeremy L Thompson 130241a4b83SYohann // Input vectors 1319e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 1322b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 1339e201c85SYohann if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 1349e201c85SYohann data->fields.inputs[i] = NULL; 135241a4b83SYohann } else { 136681d0ea7SJeremy L Thompson bool is_active; 137681d0ea7SJeremy L Thompson CeedVector vec; 138681d0ea7SJeremy L Thompson 139241a4b83SYohann // Get input vector 1402b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 141681d0ea7SJeremy L Thompson is_active = vec == CEED_VECTOR_ACTIVE; 142ea04d07fSJeremy L Thompson if (is_active) data->fields.inputs[i] = input_arr; 143ea04d07fSJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i])); 144ea04d07fSJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 145241a4b83SYohann } 146241a4b83SYohann } 147241a4b83SYohann 148241a4b83SYohann // Output vectors 1499e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 1502b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1519e201c85SYohann if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 1529e201c85SYohann data->fields.outputs[i] = NULL; 153241a4b83SYohann } else { 154681d0ea7SJeremy L Thompson bool is_active; 155681d0ea7SJeremy L Thompson CeedVector vec; 156681d0ea7SJeremy L Thompson 157241a4b83SYohann // Get output vector 1582b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 159681d0ea7SJeremy L Thompson is_active = vec == CEED_VECTOR_ACTIVE; 160ea04d07fSJeremy L Thompson if (is_active) data->fields.outputs[i] = output_arr; 1610c8fbeedSJeremy L Thompson else CeedCallBackend(CeedVectorGetArray(vec, CEED_MEM_DEVICE, &data->fields.outputs[i])); 162ea04d07fSJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 163241a4b83SYohann } 164241a4b83SYohann } 165241a4b83SYohann 1668b97b69aSJeremy L Thompson // Point coordinates, if needed 1678b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorIsAtPoints(op, &is_at_points)); 1688b97b69aSJeremy L Thompson if (is_at_points) { 1698b97b69aSJeremy L Thompson // Coords 1708b97b69aSJeremy L Thompson CeedVector vec; 1718b97b69aSJeremy L Thompson 1728b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec)); 1738b97b69aSJeremy L Thompson CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords)); 1748b97b69aSJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 1758b97b69aSJeremy L Thompson 1768b97b69aSJeremy L Thompson // Points per elem 1778b97b69aSJeremy L Thompson if (num_elem != data->points.num_elem) { 1788b97b69aSJeremy L Thompson CeedInt *points_per_elem; 1798b97b69aSJeremy L Thompson const CeedInt num_bytes = num_elem * sizeof(CeedInt); 1808b97b69aSJeremy L Thompson CeedElemRestriction rstr_points = NULL; 1818b97b69aSJeremy L Thompson 1828b97b69aSJeremy L Thompson data->points.num_elem = num_elem; 1838b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL)); 1848b97b69aSJeremy L Thompson CeedCallBackend(CeedCalloc(num_elem, &points_per_elem)); 1858b97b69aSJeremy L Thompson for (CeedInt e = 0; e < num_elem; e++) { 1868b97b69aSJeremy L Thompson CeedInt num_points_elem; 1878b97b69aSJeremy L Thompson 1888b97b69aSJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem)); 1898b97b69aSJeremy L Thompson points_per_elem[e] = num_points_elem; 1908b97b69aSJeremy L Thompson } 1918b97b69aSJeremy L Thompson if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem)); 1928b97b69aSJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes)); 1938b97b69aSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice)); 1948b97b69aSJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points)); 1958b97b69aSJeremy L Thompson CeedCallBackend(CeedFree(&points_per_elem)); 1968b97b69aSJeremy L Thompson } 1978b97b69aSJeremy L Thompson } 1988b97b69aSJeremy L Thompson 199777ff853SJeremy L Thompson // Get context data 2002b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c)); 201241a4b83SYohann 202241a4b83SYohann // Apply operator 2038b97b69aSJeremy L Thompson void *opargs[] = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W, &data->points}; 2048b97b69aSJeremy L Thompson int max_threads_per_block, min_grid_size, grid; 205ca735530SJeremy L Thompson 206dc007f05SJeremy L Thompson CeedCallBackend(CeedOperatorHasTensorBases(op, &is_tensor)); 2072b730f8bSJeremy L Thompson CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000)); 208a61b1c91SJeremy L Thompson int block[3] = {data->thread_1d, ((!is_tensor || data->dim == 1) ? 1 : data->thread_1d), -1}; 209ca735530SJeremy L Thompson 210f82027a4SJeremy L Thompson if (is_tensor) { 21190c30374SJeremy L Thompson CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, is_at_points ? 1 : max_threads_per_block, 2122b730f8bSJeremy L Thompson cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid)); 213f82027a4SJeremy L Thompson } else { 214a61b1c91SJeremy L Thompson CeedInt elems_per_block = CeedIntMin(cuda_data->device_prop.maxThreadsDim[2], CeedIntMax(512 / data->thread_1d, 1)); 215f82027a4SJeremy L Thompson 216f82027a4SJeremy L Thompson grid = num_elem / elems_per_block + (num_elem % elems_per_block > 0); 217f82027a4SJeremy L Thompson block[2] = elems_per_block; 218f82027a4SJeremy L Thompson } 21939532cebSJed Brown CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar); 220ca735530SJeremy L Thompson 221e9c76bddSJeremy L Thompson CeedCallBackend(CeedTryRunKernelDimShared_Cuda(ceed, data->op, stream, grid, block[0], block[1], block[2], shared_mem, is_run_good, opargs)); 222241a4b83SYohann 223241a4b83SYohann // Restore input arrays 2249e201c85SYohann for (CeedInt i = 0; i < num_input_fields; i++) { 2252b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 2269e201c85SYohann if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 227241a4b83SYohann } else { 228681d0ea7SJeremy L Thompson bool is_active; 229681d0ea7SJeremy L Thompson CeedVector vec; 230681d0ea7SJeremy L Thompson 2312b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 232681d0ea7SJeremy L Thompson is_active = vec == CEED_VECTOR_ACTIVE; 233ea04d07fSJeremy L Thompson if (!is_active) CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i])); 234ea04d07fSJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 235241a4b83SYohann } 236241a4b83SYohann } 237241a4b83SYohann 238241a4b83SYohann // Restore output arrays 2399e201c85SYohann for (CeedInt i = 0; i < num_output_fields; i++) { 2402b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 2419e201c85SYohann if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 242241a4b83SYohann } else { 243681d0ea7SJeremy L Thompson bool is_active; 244681d0ea7SJeremy L Thompson CeedVector vec; 245681d0ea7SJeremy L Thompson 2462b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 247681d0ea7SJeremy L Thompson is_active = vec == CEED_VECTOR_ACTIVE; 248ea04d07fSJeremy L Thompson if (!is_active) CeedCallBackend(CeedVectorRestoreArray(vec, &data->fields.outputs[i])); 249ea04d07fSJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 250241a4b83SYohann } 2513b2939feSjeremylt } 252777ff853SJeremy L Thompson 2538b97b69aSJeremy L Thompson // Restore point coordinates, if needed 2548b97b69aSJeremy L Thompson if (is_at_points) { 2558b97b69aSJeremy L Thompson CeedVector vec; 2568b97b69aSJeremy L Thompson 2578b97b69aSJeremy L Thompson CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec)); 2588b97b69aSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords)); 2598b97b69aSJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 2608b97b69aSJeremy L Thompson } 2618b97b69aSJeremy L Thompson 262777ff853SJeremy L Thompson // Restore context data 2632b730f8bSJeremy L Thompson CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c)); 264ddae5012SJeremy L Thompson 265ddae5012SJeremy L Thompson // Cleanup 2669bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 267c11e12f4SJeremy L Thompson CeedCallBackend(CeedQFunctionDestroy(&qf)); 268ea04d07fSJeremy L Thompson if (!(*is_run_good)) data->use_fallback = true; 269ea04d07fSJeremy L Thompson return CEED_ERROR_SUCCESS; 270ea04d07fSJeremy L Thompson } 271ddae5012SJeremy L Thompson 272ea04d07fSJeremy L Thompson static int CeedOperatorApplyAdd_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) { 273ea04d07fSJeremy L Thompson bool is_run_good = false; 274ea04d07fSJeremy L Thompson const CeedScalar *input_arr = NULL; 275ea04d07fSJeremy L Thompson CeedScalar *output_arr = NULL; 276ea04d07fSJeremy L Thompson 277ea04d07fSJeremy L Thompson // Try to run kernel 278ea04d07fSJeremy L Thompson if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(input_vec, CEED_MEM_DEVICE, &input_arr)); 279ea04d07fSJeremy L Thompson if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArray(output_vec, CEED_MEM_DEVICE, &output_arr)); 280e9c76bddSJeremy L Thompson CeedCallBackend(CeedOperatorApplyAddCore_Cuda_gen(op, NULL, input_arr, output_arr, &is_run_good, request)); 281ea04d07fSJeremy L Thompson if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArrayRead(input_vec, &input_arr)); 282ea04d07fSJeremy L Thompson if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArray(output_vec, &output_arr)); 283ea04d07fSJeremy L Thompson 284ea04d07fSJeremy L Thompson // Fallback on unsuccessful run 2858d12f40eSJeremy L Thompson if (!is_run_good) { 286ddae5012SJeremy L Thompson CeedOperator op_fallback; 287ddae5012SJeremy L Thompson 288ea04d07fSJeremy L Thompson CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator"); 289ddae5012SJeremy L Thompson CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback)); 290ddae5012SJeremy L Thompson CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request)); 291ddae5012SJeremy L Thompson } 292e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 293241a4b83SYohann } 294241a4b83SYohann 295c99afcd8SJeremy L Thompson static int CeedOperatorApplyAddComposite_Cuda_gen(CeedOperator op, CeedVector input_vec, CeedVector output_vec, CeedRequest *request) { 296c99afcd8SJeremy L Thompson bool is_run_good[CEED_COMPOSITE_MAX] = {false}; 297c99afcd8SJeremy L Thompson CeedInt num_suboperators; 298c99afcd8SJeremy L Thompson const CeedScalar *input_arr = NULL; 299c99afcd8SJeremy L Thompson CeedScalar *output_arr = NULL; 300087855afSJeremy L Thompson Ceed ceed; 301c99afcd8SJeremy L Thompson CeedOperator *sub_operators; 302c99afcd8SJeremy L Thompson 303087855afSJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 304c99afcd8SJeremy L Thompson CeedCall(CeedCompositeOperatorGetNumSub(op, &num_suboperators)); 305c99afcd8SJeremy L Thompson CeedCall(CeedCompositeOperatorGetSubList(op, &sub_operators)); 306c99afcd8SJeremy L Thompson if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(input_vec, CEED_MEM_DEVICE, &input_arr)); 307c99afcd8SJeremy L Thompson if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArray(output_vec, CEED_MEM_DEVICE, &output_arr)); 308c99afcd8SJeremy L Thompson for (CeedInt i = 0; i < num_suboperators; i++) { 309c99afcd8SJeremy L Thompson CeedInt num_elem = 0; 310c99afcd8SJeremy L Thompson 311c99afcd8SJeremy L Thompson CeedCall(CeedOperatorGetNumElements(sub_operators[i], &num_elem)); 312c99afcd8SJeremy L Thompson if (num_elem > 0) { 313087855afSJeremy L Thompson cudaStream_t stream = NULL; 314087855afSJeremy L Thompson 315087855afSJeremy L Thompson CeedCallCuda(ceed, cudaStreamCreate(&stream)); 316087855afSJeremy L Thompson CeedCallBackend(CeedOperatorApplyAddCore_Cuda_gen(sub_operators[i], stream, input_arr, output_arr, &is_run_good[i], request)); 317087855afSJeremy L Thompson CeedCallCuda(ceed, cudaStreamDestroy(stream)); 318c99afcd8SJeremy L Thompson } 319c99afcd8SJeremy L Thompson } 320c99afcd8SJeremy L Thompson if (input_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArrayRead(input_vec, &input_arr)); 321c99afcd8SJeremy L Thompson if (output_vec != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorRestoreArray(output_vec, &output_arr)); 322087855afSJeremy L Thompson CeedCallCuda(ceed, cudaDeviceSynchronize()); 323c99afcd8SJeremy L Thompson 324c99afcd8SJeremy L Thompson // Fallback on unsuccessful run 325c99afcd8SJeremy L Thompson for (CeedInt i = 0; i < num_suboperators; i++) { 326c99afcd8SJeremy L Thompson if (!is_run_good[i]) { 327c99afcd8SJeremy L Thompson CeedOperator op_fallback; 328c99afcd8SJeremy L Thompson 329087855afSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator"); 330c99afcd8SJeremy L Thompson CeedCallBackend(CeedOperatorGetFallback(sub_operators[i], &op_fallback)); 331c99afcd8SJeremy L Thompson CeedCallBackend(CeedOperatorApplyAdd(op_fallback, input_vec, output_vec, request)); 332c99afcd8SJeremy L Thompson } 333c99afcd8SJeremy L Thompson } 334087855afSJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 335c99afcd8SJeremy L Thompson return CEED_ERROR_SUCCESS; 336c99afcd8SJeremy L Thompson } 337c99afcd8SJeremy L Thompson 338ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 339*0183ed61SJeremy L Thompson // AtPoints diagonal assembly 340*0183ed61SJeremy L Thompson //------------------------------------------------------------------------------ 341*0183ed61SJeremy L Thompson static int CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda_gen(CeedOperator op, CeedVector assembled, CeedRequest *request) { 342*0183ed61SJeremy L Thompson Ceed ceed; 343*0183ed61SJeremy L Thompson CeedOperator_Cuda_gen *data; 344*0183ed61SJeremy L Thompson 345*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 346*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorGetData(op, &data)); 347*0183ed61SJeremy L Thompson 348*0183ed61SJeremy L Thompson // Build the assembly kernel 349*0183ed61SJeremy L Thompson if (!data->assemble_diagonal && !data->use_assembly_fallback) { 350*0183ed61SJeremy L Thompson bool is_build_good = false; 351*0183ed61SJeremy L Thompson CeedInt num_active_bases_in, num_active_bases_out; 352*0183ed61SJeremy L Thompson CeedOperatorAssemblyData assembly_data; 353*0183ed61SJeremy L Thompson 354*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorGetOperatorAssemblyData(op, &assembly_data)); 355*0183ed61SJeremy L Thompson CeedCallBackend( 356*0183ed61SJeremy L Thompson CeedOperatorAssemblyDataGetEvalModes(assembly_data, &num_active_bases_in, NULL, NULL, NULL, &num_active_bases_out, NULL, NULL, NULL, NULL)); 357*0183ed61SJeremy L Thompson if (num_active_bases_in == num_active_bases_out) { 358*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorBuildKernel_Cuda_gen(op, &is_build_good)); 359*0183ed61SJeremy L Thompson if (is_build_good) CeedCallBackend(CeedOperatorBuildKernelDiagonalAssemblyAtPoints_Cuda_gen(op, &is_build_good)); 360*0183ed61SJeremy L Thompson } 361*0183ed61SJeremy L Thompson if (!is_build_good) data->use_assembly_fallback = true; 362*0183ed61SJeremy L Thompson } 363*0183ed61SJeremy L Thompson 364*0183ed61SJeremy L Thompson // Try assembly 365*0183ed61SJeremy L Thompson if (!data->use_assembly_fallback) { 366*0183ed61SJeremy L Thompson bool is_run_good = true; 367*0183ed61SJeremy L Thompson Ceed_Cuda *cuda_data; 368*0183ed61SJeremy L Thompson CeedInt num_elem, num_input_fields, num_output_fields; 369*0183ed61SJeremy L Thompson CeedEvalMode eval_mode; 370*0183ed61SJeremy L Thompson CeedScalar *assembled_array; 371*0183ed61SJeremy L Thompson CeedQFunctionField *qf_input_fields, *qf_output_fields; 372*0183ed61SJeremy L Thompson CeedQFunction_Cuda_gen *qf_data; 373*0183ed61SJeremy L Thompson CeedQFunction qf; 374*0183ed61SJeremy L Thompson CeedOperatorField *op_input_fields, *op_output_fields; 375*0183ed61SJeremy L Thompson 376*0183ed61SJeremy L Thompson CeedCallBackend(CeedGetData(ceed, &cuda_data)); 377*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 378*0183ed61SJeremy L Thompson CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 379*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 380*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 381*0183ed61SJeremy L Thompson CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 382*0183ed61SJeremy L Thompson 383*0183ed61SJeremy L Thompson // Input vectors 384*0183ed61SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 385*0183ed61SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 386*0183ed61SJeremy L Thompson if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 387*0183ed61SJeremy L Thompson data->fields.inputs[i] = NULL; 388*0183ed61SJeremy L Thompson } else { 389*0183ed61SJeremy L Thompson bool is_active; 390*0183ed61SJeremy L Thompson CeedVector vec; 391*0183ed61SJeremy L Thompson 392*0183ed61SJeremy L Thompson // Get input vector 393*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 394*0183ed61SJeremy L Thompson is_active = vec == CEED_VECTOR_ACTIVE; 395*0183ed61SJeremy L Thompson if (is_active) data->fields.inputs[i] = NULL; 396*0183ed61SJeremy L Thompson else CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->fields.inputs[i])); 397*0183ed61SJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 398*0183ed61SJeremy L Thompson } 399*0183ed61SJeremy L Thompson } 400*0183ed61SJeremy L Thompson 401*0183ed61SJeremy L Thompson // Point coordinates 402*0183ed61SJeremy L Thompson { 403*0183ed61SJeremy L Thompson CeedVector vec; 404*0183ed61SJeremy L Thompson 405*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec)); 406*0183ed61SJeremy L Thompson CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &data->points.coords)); 407*0183ed61SJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 408*0183ed61SJeremy L Thompson 409*0183ed61SJeremy L Thompson // Points per elem 410*0183ed61SJeremy L Thompson if (num_elem != data->points.num_elem) { 411*0183ed61SJeremy L Thompson CeedInt *points_per_elem; 412*0183ed61SJeremy L Thompson const CeedInt num_bytes = num_elem * sizeof(CeedInt); 413*0183ed61SJeremy L Thompson CeedElemRestriction rstr_points = NULL; 414*0183ed61SJeremy L Thompson 415*0183ed61SJeremy L Thompson data->points.num_elem = num_elem; 416*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL)); 417*0183ed61SJeremy L Thompson CeedCallBackend(CeedCalloc(num_elem, &points_per_elem)); 418*0183ed61SJeremy L Thompson for (CeedInt e = 0; e < num_elem; e++) { 419*0183ed61SJeremy L Thompson CeedInt num_points_elem; 420*0183ed61SJeremy L Thompson 421*0183ed61SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem)); 422*0183ed61SJeremy L Thompson points_per_elem[e] = num_points_elem; 423*0183ed61SJeremy L Thompson } 424*0183ed61SJeremy L Thompson if (data->points.num_per_elem) CeedCallCuda(ceed, cudaFree((void **)data->points.num_per_elem)); 425*0183ed61SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&data->points.num_per_elem, num_bytes)); 426*0183ed61SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((void *)data->points.num_per_elem, points_per_elem, num_bytes, cudaMemcpyHostToDevice)); 427*0183ed61SJeremy L Thompson CeedCallBackend(CeedElemRestrictionDestroy(&rstr_points)); 428*0183ed61SJeremy L Thompson CeedCallBackend(CeedFree(&points_per_elem)); 429*0183ed61SJeremy L Thompson } 430*0183ed61SJeremy L Thompson } 431*0183ed61SJeremy L Thompson 432*0183ed61SJeremy L Thompson // Get context data 433*0183ed61SJeremy L Thompson CeedCallBackend(CeedQFunctionGetInnerContextData(qf, CEED_MEM_DEVICE, &qf_data->d_c)); 434*0183ed61SJeremy L Thompson 435*0183ed61SJeremy L Thompson // Assembly array 436*0183ed61SJeremy L Thompson CeedCallBackend(CeedVectorGetArray(assembled, CEED_MEM_DEVICE, &assembled_array)); 437*0183ed61SJeremy L Thompson 438*0183ed61SJeremy L Thompson // Assemble diagonal 439*0183ed61SJeremy L Thompson void *opargs[] = {(void *)&num_elem, &qf_data->d_c, &data->indices, &data->fields, &data->B, &data->G, &data->W, &data->points, &assembled_array}; 440*0183ed61SJeremy L Thompson int max_threads_per_block, min_grid_size, grid; 441*0183ed61SJeremy L Thompson 442*0183ed61SJeremy L Thompson CeedCallCuda(ceed, cuOccupancyMaxPotentialBlockSize(&min_grid_size, &max_threads_per_block, data->op, dynamicSMemSize, 0, 0x10000)); 443*0183ed61SJeremy L Thompson int block[3] = {data->thread_1d, (data->dim == 1 ? 1 : data->thread_1d), -1}; 444*0183ed61SJeremy L Thompson 445*0183ed61SJeremy L Thompson CeedCallBackend(BlockGridCalculate(num_elem, min_grid_size / cuda_data->device_prop.multiProcessorCount, 1, 446*0183ed61SJeremy L Thompson cuda_data->device_prop.maxThreadsDim[2], cuda_data->device_prop.warpSize, block, &grid)); 447*0183ed61SJeremy L Thompson CeedInt shared_mem = block[0] * block[1] * block[2] * sizeof(CeedScalar); 448*0183ed61SJeremy L Thompson 449*0183ed61SJeremy L Thompson CeedCallBackend( 450*0183ed61SJeremy L Thompson CeedTryRunKernelDimShared_Cuda(ceed, data->assemble_diagonal, NULL, grid, block[0], block[1], block[2], shared_mem, &is_run_good, opargs)); 451*0183ed61SJeremy L Thompson 452*0183ed61SJeremy L Thompson // Restore input arrays 453*0183ed61SJeremy L Thompson for (CeedInt i = 0; i < num_input_fields; i++) { 454*0183ed61SJeremy L Thompson CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 455*0183ed61SJeremy L Thompson if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 456*0183ed61SJeremy L Thompson } else { 457*0183ed61SJeremy L Thompson bool is_active; 458*0183ed61SJeremy L Thompson CeedVector vec; 459*0183ed61SJeremy L Thompson 460*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 461*0183ed61SJeremy L Thompson is_active = vec == CEED_VECTOR_ACTIVE; 462*0183ed61SJeremy L Thompson if (!is_active) CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->fields.inputs[i])); 463*0183ed61SJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 464*0183ed61SJeremy L Thompson } 465*0183ed61SJeremy L Thompson } 466*0183ed61SJeremy L Thompson 467*0183ed61SJeremy L Thompson // Restore point coordinates 468*0183ed61SJeremy L Thompson { 469*0183ed61SJeremy L Thompson CeedVector vec; 470*0183ed61SJeremy L Thompson 471*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorAtPointsGetPoints(op, NULL, &vec)); 472*0183ed61SJeremy L Thompson CeedCallBackend(CeedVectorRestoreArrayRead(vec, &data->points.coords)); 473*0183ed61SJeremy L Thompson CeedCallBackend(CeedVectorDestroy(&vec)); 474*0183ed61SJeremy L Thompson } 475*0183ed61SJeremy L Thompson 476*0183ed61SJeremy L Thompson // Restore context data 477*0183ed61SJeremy L Thompson CeedCallBackend(CeedQFunctionRestoreInnerContextData(qf, &qf_data->d_c)); 478*0183ed61SJeremy L Thompson 479*0183ed61SJeremy L Thompson // Restore assembly array 480*0183ed61SJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(assembled, &assembled_array)); 481*0183ed61SJeremy L Thompson 482*0183ed61SJeremy L Thompson // Cleanup 483*0183ed61SJeremy L Thompson CeedCallBackend(CeedQFunctionDestroy(&qf)); 484*0183ed61SJeremy L Thompson if (!is_run_good) data->use_assembly_fallback = true; 485*0183ed61SJeremy L Thompson } 486*0183ed61SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 487*0183ed61SJeremy L Thompson 488*0183ed61SJeremy L Thompson // Fallback, if needed 489*0183ed61SJeremy L Thompson if (data->use_assembly_fallback) { 490*0183ed61SJeremy L Thompson CeedOperator op_fallback; 491*0183ed61SJeremy L Thompson 492*0183ed61SJeremy L Thompson CeedDebug256(CeedOperatorReturnCeed(op), CEED_DEBUG_COLOR_SUCCESS, "Falling back to /gpu/cuda/ref CeedOperator"); 493*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorGetFallback(op, &op_fallback)); 494*0183ed61SJeremy L Thompson CeedCallBackend(CeedOperatorLinearAssembleAddDiagonal(op_fallback, assembled, request)); 495*0183ed61SJeremy L Thompson return CEED_ERROR_SUCCESS; 496*0183ed61SJeremy L Thompson } 497*0183ed61SJeremy L Thompson return CEED_ERROR_SUCCESS; 498*0183ed61SJeremy L Thompson } 499*0183ed61SJeremy L Thompson 500*0183ed61SJeremy L Thompson //------------------------------------------------------------------------------ 501ab213215SJeremy L Thompson // Create operator 502ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 503241a4b83SYohann int CeedOperatorCreate_Cuda_gen(CeedOperator op) { 504*0183ed61SJeremy L Thompson bool is_composite, is_at_points; 505241a4b83SYohann Ceed ceed; 506241a4b83SYohann CeedOperator_Cuda_gen *impl; 507241a4b83SYohann 508ca735530SJeremy L Thompson CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 5092b730f8bSJeremy L Thompson CeedCallBackend(CeedCalloc(1, &impl)); 5102b730f8bSJeremy L Thompson CeedCallBackend(CeedOperatorSetData(op, impl)); 511c99afcd8SJeremy L Thompson CeedCall(CeedOperatorIsComposite(op, &is_composite)); 512c99afcd8SJeremy L Thompson if (is_composite) { 513c99afcd8SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAddComposite", CeedOperatorApplyAddComposite_Cuda_gen)); 514c99afcd8SJeremy L Thompson } else { 5152b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda_gen)); 516c99afcd8SJeremy L Thompson } 517*0183ed61SJeremy L Thompson CeedCall(CeedOperatorIsAtPoints(op, &is_at_points)); 518*0183ed61SJeremy L Thompson if (is_at_points) { 519*0183ed61SJeremy L Thompson CeedCallBackend( 520*0183ed61SJeremy L Thompson CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda_gen)); 521*0183ed61SJeremy L Thompson } 5222b730f8bSJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda_gen)); 5239bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 524e15f9bd0SJeremy L Thompson return CEED_ERROR_SUCCESS; 525241a4b83SYohann } 5266aa95790SJeremy L Thompson 527ab213215SJeremy L Thompson //------------------------------------------------------------------------------ 528