1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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. 30d0321e0SJeremy L Thompson // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 50d0321e0SJeremy L Thompson // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 70d0321e0SJeremy L Thompson 849aac155SJeremy L Thompson #include <ceed.h> 90d0321e0SJeremy L Thompson #include <hip/hip_runtime.h> 100d0321e0SJeremy L Thompson 110d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 123196072fSJeremy L Thompson // Kernel for copy strided on device 133196072fSJeremy L Thompson //------------------------------------------------------------------------------ 143196072fSJeremy L Thompson __global__ static void copyStridedK(CeedScalar *__restrict__ vec, CeedSize start, CeedSize step, CeedSize size, CeedScalar *__restrict__ vec_copy) { 15c9d5affaSJeremy L Thompson const CeedSize index = threadIdx.x + (CeedSize)blockDim.x * blockIdx.x; 163196072fSJeremy L Thompson 17c9d5affaSJeremy L Thompson if (index < size) { 183196072fSJeremy L Thompson if ((index - start) % step == 0) vec_copy[index] = vec[index]; 193196072fSJeremy L Thompson } 20c9d5affaSJeremy L Thompson } 213196072fSJeremy L Thompson 223196072fSJeremy L Thompson //------------------------------------------------------------------------------ 233196072fSJeremy L Thompson // Copy strided on device memory 243196072fSJeremy L Thompson //------------------------------------------------------------------------------ 253196072fSJeremy L Thompson extern "C" int CeedDeviceCopyStrided_Hip(CeedScalar *d_array, CeedSize start, CeedSize step, CeedSize length, CeedScalar *d_copy_array) { 263196072fSJeremy L Thompson const int block_size = 512; 273196072fSJeremy L Thompson const CeedSize vec_size = length; 283196072fSJeremy L Thompson int grid_size = vec_size / block_size; 293196072fSJeremy L Thompson 303196072fSJeremy L Thompson if (block_size * grid_size < vec_size) grid_size += 1; 313196072fSJeremy L Thompson hipLaunchKernelGGL(copyStridedK, dim3(grid_size), dim3(block_size), 0, 0, d_array, start, step, length, d_copy_array); 323196072fSJeremy L Thompson return 0; 333196072fSJeremy L Thompson } 343196072fSJeremy L Thompson 353196072fSJeremy L Thompson //------------------------------------------------------------------------------ 360d0321e0SJeremy L Thompson // Kernel for set value on device 370d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 389330daecSnbeams __global__ static void setValueK(CeedScalar *__restrict__ vec, CeedSize size, CeedScalar val) { 39c9d5affaSJeremy L Thompson const CeedSize index = threadIdx.x + (CeedSize)blockDim.x * blockIdx.x; 40b7453713SJeremy L Thompson 41c9d5affaSJeremy L Thompson if (index < size) vec[index] = val; 420d0321e0SJeremy L Thompson } 430d0321e0SJeremy L Thompson 440d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 450d0321e0SJeremy L Thompson // Set value on device memory 460d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 479330daecSnbeams extern "C" int CeedDeviceSetValue_Hip(CeedScalar *d_array, CeedSize length, CeedScalar val) { 48b7453713SJeremy L Thompson const int block_size = 512; 49b7453713SJeremy L Thompson const CeedSize vec_size = length; 50b7453713SJeremy L Thompson int grid_size = vec_size / block_size; 510d0321e0SJeremy L Thompson 52b7453713SJeremy L Thompson if (block_size * grid_size < vec_size) grid_size += 1; 53b7453713SJeremy L Thompson hipLaunchKernelGGL(setValueK, dim3(grid_size), dim3(block_size), 0, 0, d_array, length, val); 540d0321e0SJeremy L Thompson return 0; 550d0321e0SJeremy L Thompson } 560d0321e0SJeremy L Thompson 570d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 583196072fSJeremy L Thompson // Kernel for set value strided on device 593196072fSJeremy L Thompson //------------------------------------------------------------------------------ 602d73a370SJeremy L Thompson __global__ static void setValueStridedK(CeedScalar *__restrict__ vec, CeedSize start, CeedSize stop, CeedSize step, CeedScalar val) { 61c9d5affaSJeremy L Thompson const CeedSize index = threadIdx.x + (CeedSize)blockDim.x * blockIdx.x; 623196072fSJeremy L Thompson 6314c82621SJeremy L Thompson if (index < stop - start) { 6414c82621SJeremy L Thompson if (index % step == 0) vec[start + index] = val; 653196072fSJeremy L Thompson } 66c9d5affaSJeremy L Thompson } 673196072fSJeremy L Thompson 683196072fSJeremy L Thompson //------------------------------------------------------------------------------ 693196072fSJeremy L Thompson // Set value strided on device memory 703196072fSJeremy L Thompson //------------------------------------------------------------------------------ 71ff90b007SJeremy L Thompson extern "C" int CeedDeviceSetValueStrided_Hip(CeedScalar *d_array, CeedSize start, CeedInt stop, CeedSize step, CeedSize length, CeedScalar val) { 723196072fSJeremy L Thompson const int block_size = 512; 7314c82621SJeremy L Thompson const CeedSize set_size = stop - start; 7414c82621SJeremy L Thompson int grid_size = set_size / block_size; 753196072fSJeremy L Thompson 7614c82621SJeremy L Thompson if (block_size * grid_size < set_size) grid_size += 1; 772d73a370SJeremy L Thompson hipLaunchKernelGGL(setValueStridedK, dim3(grid_size), dim3(block_size), 0, 0, d_array, start, stop, step, val); 783196072fSJeremy L Thompson return 0; 793196072fSJeremy L Thompson } 803196072fSJeremy L Thompson 813196072fSJeremy L Thompson //------------------------------------------------------------------------------ 820d0321e0SJeremy L Thompson // Kernel for taking reciprocal 830d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 849330daecSnbeams __global__ static void rcpValueK(CeedScalar *__restrict__ vec, CeedSize size) { 85c9d5affaSJeremy L Thompson const CeedSize index = threadIdx.x + (CeedSize)blockDim.x * blockIdx.x; 86b7453713SJeremy L Thompson 87c9d5affaSJeremy L Thompson if (index < size) { 88b7453713SJeremy L Thompson if (fabs(vec[index]) > 1E-16) vec[index] = 1. / vec[index]; 890d0321e0SJeremy L Thompson } 90c9d5affaSJeremy L Thompson } 910d0321e0SJeremy L Thompson 920d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 930d0321e0SJeremy L Thompson // Take vector reciprocal in device memory 940d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 959330daecSnbeams extern "C" int CeedDeviceReciprocal_Hip(CeedScalar *d_array, CeedSize length) { 96b7453713SJeremy L Thompson const int block_size = 512; 97b7453713SJeremy L Thompson const CeedSize vec_size = length; 98b7453713SJeremy L Thompson int grid_size = vec_size / block_size; 990d0321e0SJeremy L Thompson 100b7453713SJeremy L Thompson if (block_size * grid_size < vec_size) grid_size += 1; 101b7453713SJeremy L Thompson hipLaunchKernelGGL(rcpValueK, dim3(grid_size), dim3(block_size), 0, 0, d_array, length); 1020d0321e0SJeremy L Thompson return 0; 1030d0321e0SJeremy L Thompson } 1040d0321e0SJeremy L Thompson 1050d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1060d0321e0SJeremy L Thompson // Kernel for scale 1070d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1089330daecSnbeams __global__ static void scaleValueK(CeedScalar *__restrict__ x, CeedScalar alpha, CeedSize size) { 109c9d5affaSJeremy L Thompson const CeedSize index = threadIdx.x + (CeedSize)blockDim.x * blockIdx.x; 110b7453713SJeremy L Thompson 111c9d5affaSJeremy L Thompson if (index < size) x[index] *= alpha; 1120d0321e0SJeremy L Thompson } 1130d0321e0SJeremy L Thompson 1140d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1150d0321e0SJeremy L Thompson // Compute x = alpha x on device 1160d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1179330daecSnbeams extern "C" int CeedDeviceScale_Hip(CeedScalar *x_array, CeedScalar alpha, CeedSize length) { 118b7453713SJeremy L Thompson const int block_size = 512; 119b7453713SJeremy L Thompson const CeedSize vec_size = length; 120b7453713SJeremy L Thompson int grid_size = vec_size / block_size; 1210d0321e0SJeremy L Thompson 122b7453713SJeremy L Thompson if (block_size * grid_size < vec_size) grid_size += 1; 123b7453713SJeremy L Thompson hipLaunchKernelGGL(scaleValueK, dim3(grid_size), dim3(block_size), 0, 0, x_array, alpha, length); 1240d0321e0SJeremy L Thompson return 0; 1250d0321e0SJeremy L Thompson } 1260d0321e0SJeremy L Thompson 1270d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1280d0321e0SJeremy L Thompson // Kernel for axpy 1290d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1309330daecSnbeams __global__ static void axpyValueK(CeedScalar *__restrict__ y, CeedScalar alpha, CeedScalar *__restrict__ x, CeedSize size) { 131c9d5affaSJeremy L Thompson const CeedSize index = threadIdx.x + (CeedSize)blockDim.x * blockIdx.x; 132b7453713SJeremy L Thompson 133c9d5affaSJeremy L Thompson if (index < size) y[index] += alpha * x[index]; 1340d0321e0SJeremy L Thompson } 1350d0321e0SJeremy L Thompson 1360d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1370d0321e0SJeremy L Thompson // Compute y = alpha x + y on device 1380d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1399330daecSnbeams extern "C" int CeedDeviceAXPY_Hip(CeedScalar *y_array, CeedScalar alpha, CeedScalar *x_array, CeedSize length) { 140b7453713SJeremy L Thompson const int block_size = 512; 141b7453713SJeremy L Thompson const CeedSize vec_size = length; 142b7453713SJeremy L Thompson int grid_size = vec_size / block_size; 1430d0321e0SJeremy L Thompson 144b7453713SJeremy L Thompson if (block_size * grid_size < vec_size) grid_size += 1; 145b7453713SJeremy L Thompson hipLaunchKernelGGL(axpyValueK, dim3(grid_size), dim3(block_size), 0, 0, y_array, alpha, x_array, length); 1460d0321e0SJeremy L Thompson return 0; 1470d0321e0SJeremy L Thompson } 1480d0321e0SJeremy L Thompson 1490d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1505fb68f37SKaren (Ren) Stengel // Kernel for axpby 1515fb68f37SKaren (Ren) Stengel //------------------------------------------------------------------------------ 1529330daecSnbeams __global__ static void axpbyValueK(CeedScalar *__restrict__ y, CeedScalar alpha, CeedScalar beta, CeedScalar *__restrict__ x, CeedSize size) { 153c9d5affaSJeremy L Thompson const CeedSize index = threadIdx.x + (CeedSize)blockDim.x * blockIdx.x; 154b7453713SJeremy L Thompson 155c9d5affaSJeremy L Thompson if (index < size) { 156b7453713SJeremy L Thompson y[index] = beta * y[index]; 157b7453713SJeremy L Thompson y[index] += alpha * x[index]; 1585fb68f37SKaren (Ren) Stengel } 159c9d5affaSJeremy L Thompson } 1605fb68f37SKaren (Ren) Stengel 1615fb68f37SKaren (Ren) Stengel //------------------------------------------------------------------------------ 1625fb68f37SKaren (Ren) Stengel // Compute y = alpha x + beta y on device 1635fb68f37SKaren (Ren) Stengel //------------------------------------------------------------------------------ 1649330daecSnbeams extern "C" int CeedDeviceAXPBY_Hip(CeedScalar *y_array, CeedScalar alpha, CeedScalar beta, CeedScalar *x_array, CeedSize length) { 165b7453713SJeremy L Thompson const int block_size = 512; 166b7453713SJeremy L Thompson const CeedSize vec_size = length; 167b7453713SJeremy L Thompson int grid_size = vec_size / block_size; 1685fb68f37SKaren (Ren) Stengel 169b7453713SJeremy L Thompson if (block_size * grid_size < vec_size) grid_size += 1; 170b7453713SJeremy L Thompson hipLaunchKernelGGL(axpbyValueK, dim3(grid_size), dim3(block_size), 0, 0, y_array, alpha, beta, x_array, length); 1715fb68f37SKaren (Ren) Stengel return 0; 1725fb68f37SKaren (Ren) Stengel } 1735fb68f37SKaren (Ren) Stengel 1745fb68f37SKaren (Ren) Stengel //------------------------------------------------------------------------------ 1750d0321e0SJeremy L Thompson // Kernel for pointwise mult 1760d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1779330daecSnbeams __global__ static void pointwiseMultValueK(CeedScalar *__restrict__ w, CeedScalar *x, CeedScalar *__restrict__ y, CeedSize size) { 178c9d5affaSJeremy L Thompson const CeedSize index = threadIdx.x + (CeedSize)blockDim.x * blockIdx.x; 179b7453713SJeremy L Thompson 180c9d5affaSJeremy L Thompson if (index < size) w[index] = x[index] * y[index]; 1810d0321e0SJeremy L Thompson } 1820d0321e0SJeremy L Thompson 1830d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1840d0321e0SJeremy L Thompson // Compute the pointwise multiplication w = x .* y on device 1850d0321e0SJeremy L Thompson //------------------------------------------------------------------------------ 1869330daecSnbeams extern "C" int CeedDevicePointwiseMult_Hip(CeedScalar *w_array, CeedScalar *x_array, CeedScalar *y_array, CeedSize length) { 187b7453713SJeremy L Thompson const int block_size = 512; 188b7453713SJeremy L Thompson const CeedSize vec_size = length; 189b7453713SJeremy L Thompson int grid_size = vec_size / block_size; 1900d0321e0SJeremy L Thompson 191b7453713SJeremy L Thompson if (block_size * grid_size < vec_size) grid_size += 1; 192b7453713SJeremy L Thompson hipLaunchKernelGGL(pointwiseMultValueK, dim3(grid_size), dim3(block_size), 0, 0, w_array, x_array, y_array, length); 1930d0321e0SJeremy L Thompson return 0; 1940d0321e0SJeremy L Thompson } 1952a86cc9dSSebastian Grimberg 1962a86cc9dSSebastian Grimberg //------------------------------------------------------------------------------ 197