xref: /libCEED/backends/cuda-ref/kernels/cuda-ref-vector.cu (revision 2a86cc9d4dbfce2964c7e8927a1e6db8d19a41fc)
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.
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 <cuda.h>
100d0321e0SJeremy L Thompson 
110d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
120d0321e0SJeremy L Thompson // Kernel for set value on device
130d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
140d0321e0SJeremy L Thompson __global__ static void setValueK(CeedScalar * __restrict__ vec, CeedInt size,
150d0321e0SJeremy L Thompson                                  CeedScalar val) {
160d0321e0SJeremy L Thompson   int idx = threadIdx.x + blockDim.x * blockIdx.x;
170d0321e0SJeremy L Thompson   if (idx >= size)
180d0321e0SJeremy L Thompson     return;
190d0321e0SJeremy L Thompson   vec[idx] = val;
200d0321e0SJeremy L Thompson }
210d0321e0SJeremy L Thompson 
220d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
230d0321e0SJeremy L Thompson // Set value on device memory
240d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
250d0321e0SJeremy L Thompson extern "C" int CeedDeviceSetValue_Cuda(CeedScalar* d_array, CeedInt length,
260d0321e0SJeremy L Thompson                                        CeedScalar val) {
270d0321e0SJeremy L Thompson   const int bsize = 512;
280d0321e0SJeremy L Thompson   const int vecsize = length;
290d0321e0SJeremy L Thompson   int gridsize = vecsize / bsize;
300d0321e0SJeremy L Thompson 
310d0321e0SJeremy L Thompson   if (bsize * gridsize < vecsize)
320d0321e0SJeremy L Thompson     gridsize += 1;
330d0321e0SJeremy L Thompson   setValueK<<<gridsize,bsize>>>(d_array, length, val);
340d0321e0SJeremy L Thompson   return 0;
350d0321e0SJeremy L Thompson }
360d0321e0SJeremy L Thompson 
370d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
380d0321e0SJeremy L Thompson // Kernel for taking reciprocal
390d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
400d0321e0SJeremy L Thompson __global__ static void rcpValueK(CeedScalar * __restrict__ vec, CeedInt size) {
410d0321e0SJeremy L Thompson   int idx = threadIdx.x + blockDim.x * blockIdx.x;
420d0321e0SJeremy L Thompson   if (idx >= size)
430d0321e0SJeremy L Thompson     return;
440d0321e0SJeremy L Thompson   if (fabs(vec[idx]) > 1E-16)
450d0321e0SJeremy L Thompson     vec[idx] = 1./vec[idx];
460d0321e0SJeremy L Thompson }
470d0321e0SJeremy L Thompson 
480d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
490d0321e0SJeremy L Thompson // Take vector reciprocal in device memory
500d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
510d0321e0SJeremy L Thompson extern "C" int CeedDeviceReciprocal_Cuda(CeedScalar* d_array, CeedInt length) {
520d0321e0SJeremy L Thompson   const int bsize = 512;
530d0321e0SJeremy L Thompson   const int vecsize = length;
540d0321e0SJeremy L Thompson   int gridsize = vecsize / bsize;
550d0321e0SJeremy L Thompson 
560d0321e0SJeremy L Thompson   if (bsize * gridsize < vecsize)
570d0321e0SJeremy L Thompson     gridsize += 1;
580d0321e0SJeremy L Thompson   rcpValueK<<<gridsize,bsize>>>(d_array, length);
590d0321e0SJeremy L Thompson   return 0;
600d0321e0SJeremy L Thompson }
610d0321e0SJeremy L Thompson 
620d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
630d0321e0SJeremy L Thompson // Kernel for scale
640d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
650d0321e0SJeremy L Thompson __global__ static void scaleValueK(CeedScalar * __restrict__ x, CeedScalar alpha,
660d0321e0SJeremy L Thompson     CeedInt size) {
670d0321e0SJeremy L Thompson   int idx = threadIdx.x + blockDim.x * blockIdx.x;
680d0321e0SJeremy L Thompson   if (idx >= size)
690d0321e0SJeremy L Thompson     return;
700d0321e0SJeremy L Thompson   x[idx] *= alpha;
710d0321e0SJeremy L Thompson }
720d0321e0SJeremy L Thompson 
730d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
740d0321e0SJeremy L Thompson // Compute x = alpha x on device
750d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
760d0321e0SJeremy L Thompson extern "C" int CeedDeviceScale_Cuda(CeedScalar *x_array, CeedScalar alpha,
770d0321e0SJeremy L Thompson     CeedInt length) {
780d0321e0SJeremy L Thompson   const int bsize = 512;
790d0321e0SJeremy L Thompson   const int vecsize = length;
800d0321e0SJeremy L Thompson   int gridsize = vecsize / bsize;
810d0321e0SJeremy L Thompson 
820d0321e0SJeremy L Thompson   if (bsize * gridsize < vecsize)
830d0321e0SJeremy L Thompson     gridsize += 1;
840d0321e0SJeremy L Thompson   scaleValueK<<<gridsize,bsize>>>(x_array, alpha, length);
850d0321e0SJeremy L Thompson   return 0;
860d0321e0SJeremy L Thompson }
870d0321e0SJeremy L Thompson 
880d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
890d0321e0SJeremy L Thompson // Kernel for axpy
900d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
910d0321e0SJeremy L Thompson __global__ static void axpyValueK(CeedScalar * __restrict__ y, CeedScalar alpha,
920d0321e0SJeremy L Thompson     CeedScalar * __restrict__ x, CeedInt size) {
930d0321e0SJeremy L Thompson   int idx = threadIdx.x + blockDim.x * blockIdx.x;
940d0321e0SJeremy L Thompson   if (idx >= size)
950d0321e0SJeremy L Thompson     return;
960d0321e0SJeremy L Thompson   y[idx] += alpha * x[idx];
970d0321e0SJeremy L Thompson }
980d0321e0SJeremy L Thompson 
990d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1000d0321e0SJeremy L Thompson // Compute y = alpha x + y on device
1010d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1020d0321e0SJeremy L Thompson extern "C" int CeedDeviceAXPY_Cuda(CeedScalar *y_array, CeedScalar alpha,
1030d0321e0SJeremy L Thompson     CeedScalar *x_array, CeedInt length) {
1040d0321e0SJeremy L Thompson   const int bsize = 512;
1050d0321e0SJeremy L Thompson   const int vecsize = length;
1060d0321e0SJeremy L Thompson   int gridsize = vecsize / bsize;
1070d0321e0SJeremy L Thompson 
1080d0321e0SJeremy L Thompson   if (bsize * gridsize < vecsize)
1090d0321e0SJeremy L Thompson     gridsize += 1;
1100d0321e0SJeremy L Thompson   axpyValueK<<<gridsize,bsize>>>(y_array, alpha, x_array, length);
1110d0321e0SJeremy L Thompson   return 0;
1120d0321e0SJeremy L Thompson }
1130d0321e0SJeremy L Thompson 
1140d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1155fb68f37SKaren (Ren) Stengel // Kernel for axpby
1165fb68f37SKaren (Ren) Stengel //------------------------------------------------------------------------------
1175fb68f37SKaren (Ren) Stengel __global__ static void axpbyValueK(CeedScalar * __restrict__ y, CeedScalar alpha, CeedScalar beta,
1185fb68f37SKaren (Ren) Stengel     CeedScalar * __restrict__ x, CeedInt size) {
1195fb68f37SKaren (Ren) Stengel   int idx = threadIdx.x + blockDim.x * blockIdx.x;
1205fb68f37SKaren (Ren) Stengel   if (idx >= size)
1215fb68f37SKaren (Ren) Stengel     return;
1225fb68f37SKaren (Ren) Stengel   y[idx] = beta * y[idx];
1235fb68f37SKaren (Ren) Stengel   y[idx] += alpha * x[idx];
1245fb68f37SKaren (Ren) Stengel }
1255fb68f37SKaren (Ren) Stengel 
1265fb68f37SKaren (Ren) Stengel //------------------------------------------------------------------------------
1275fb68f37SKaren (Ren) Stengel // Compute y = alpha x + beta y on device
1285fb68f37SKaren (Ren) Stengel //------------------------------------------------------------------------------
1295fb68f37SKaren (Ren) Stengel extern "C" int CeedDeviceAXPBY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar beta,
1305fb68f37SKaren (Ren) Stengel     CeedScalar *x_array, CeedInt length) {
1315fb68f37SKaren (Ren) Stengel   const int bsize = 512;
1325fb68f37SKaren (Ren) Stengel   const int vecsize = length;
1335fb68f37SKaren (Ren) Stengel   int gridsize = vecsize / bsize;
1345fb68f37SKaren (Ren) Stengel 
1355fb68f37SKaren (Ren) Stengel   if (bsize * gridsize < vecsize)
1365fb68f37SKaren (Ren) Stengel     gridsize += 1;
1375fb68f37SKaren (Ren) Stengel   axpbyValueK<<<gridsize,bsize>>>(y_array, alpha, beta, x_array, length);
1385fb68f37SKaren (Ren) Stengel   return 0;
1395fb68f37SKaren (Ren) Stengel }
1405fb68f37SKaren (Ren) Stengel 
1415fb68f37SKaren (Ren) Stengel //------------------------------------------------------------------------------
1420d0321e0SJeremy L Thompson // Kernel for pointwise mult
1430d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1440d0321e0SJeremy L Thompson __global__ static void pointwiseMultValueK(CeedScalar * __restrict__ w,
1450d0321e0SJeremy L Thompson     CeedScalar * x, CeedScalar * __restrict__ y, CeedInt size) {
1460d0321e0SJeremy L Thompson   int idx = threadIdx.x + blockDim.x * blockIdx.x;
1470d0321e0SJeremy L Thompson   if (idx >= size)
1480d0321e0SJeremy L Thompson     return;
1490d0321e0SJeremy L Thompson   w[idx] = x[idx] * y[idx];
1500d0321e0SJeremy L Thompson }
1510d0321e0SJeremy L Thompson 
1520d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1530d0321e0SJeremy L Thompson // Compute the pointwise multiplication w = x .* y on device
1540d0321e0SJeremy L Thompson //------------------------------------------------------------------------------
1550d0321e0SJeremy L Thompson extern "C" int CeedDevicePointwiseMult_Cuda(CeedScalar *w_array, CeedScalar *x_array,
1560d0321e0SJeremy L Thompson     CeedScalar *y_array, CeedInt length) {
1570d0321e0SJeremy L Thompson   const int bsize = 512;
1580d0321e0SJeremy L Thompson   const int vecsize = length;
1590d0321e0SJeremy L Thompson   int gridsize = vecsize / bsize;
1600d0321e0SJeremy L Thompson 
1610d0321e0SJeremy L Thompson   if (bsize * gridsize < vecsize)
1620d0321e0SJeremy L Thompson     gridsize += 1;
1630d0321e0SJeremy L Thompson   pointwiseMultValueK<<<gridsize,bsize>>>(w_array, x_array, y_array, length);
1640d0321e0SJeremy L Thompson   return 0;
1650d0321e0SJeremy L Thompson }
166*2a86cc9dSSebastian Grimberg 
167*2a86cc9dSSebastian Grimberg //------------------------------------------------------------------------------
168