xref: /petsc/src/mat/impls/sell/seq/seqcuda/sellcuda.cu (revision 8df136f97abf49e837a1a7f058b7329c31cd7d87)
12d1451d4SHong Zhang #include <cuda_runtime.h>
22d1451d4SHong Zhang 
32d1451d4SHong Zhang #include <petscdevice_cuda.h>
42d1451d4SHong Zhang #include <../src/mat/impls/sell/seq/sell.h> /*I   "petscmat.h"  I*/
52d1451d4SHong Zhang 
607e43b41SHong Zhang #define SLICE_HEIGHT 16
707e43b41SHong Zhang 
82d1451d4SHong Zhang typedef struct {
990d2215bSHong Zhang   PetscInt   maxallocmat;
1090d2215bSHong Zhang   PetscInt   totalentries;
1190d2215bSHong Zhang   PetscInt  *colidx; /* column index array, device pointer */
1290d2215bSHong Zhang   MatScalar *val;    /* value array, device pointer */
1390d2215bSHong Zhang   PetscInt   totalslices;
1490d2215bSHong Zhang   PetscInt  *sliidx; /* slice index array, device pointer */
152d1451d4SHong Zhang   PetscInt   nonzerostate;
1607e43b41SHong Zhang   PetscInt   kernelchoice;
174e58db63SHong Zhang   PetscInt   blocky;
1890d2215bSHong Zhang   PetscInt   chunksperblock;
1990d2215bSHong Zhang   PetscInt   totalchunks;
2090d2215bSHong Zhang   PetscInt  *chunk_slice_map; /* starting slice for each chunk, device pointer */
212d1451d4SHong Zhang } Mat_SeqSELLCUDA;
222d1451d4SHong Zhang 
232d1451d4SHong Zhang static PetscErrorCode MatSeqSELLCUDA_Destroy(Mat_SeqSELLCUDA **cudastruct)
242d1451d4SHong Zhang {
252d1451d4SHong Zhang   PetscFunctionBegin;
262d1451d4SHong Zhang   if (*cudastruct) {
272d1451d4SHong Zhang     if ((*cudastruct)->colidx) { PetscCallCUDA(cudaFree((*cudastruct)->colidx)); }
282d1451d4SHong Zhang     if ((*cudastruct)->val) { PetscCallCUDA(cudaFree((*cudastruct)->val)); }
292d1451d4SHong Zhang     if ((*cudastruct)->sliidx) { PetscCallCUDA(cudaFree((*cudastruct)->sliidx)); }
3090d2215bSHong Zhang     if ((*cudastruct)->chunk_slice_map) { PetscCallCUDA(cudaFree((*cudastruct)->chunk_slice_map)); }
312d1451d4SHong Zhang     PetscCall(PetscFree(*cudastruct));
322d1451d4SHong Zhang   }
332d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
342d1451d4SHong Zhang }
352d1451d4SHong Zhang 
362d1451d4SHong Zhang static PetscErrorCode MatSeqSELLCUDACopyToGPU(Mat A)
372d1451d4SHong Zhang {
382d1451d4SHong Zhang   Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr;
392d1451d4SHong Zhang   Mat_SeqSELL     *a          = (Mat_SeqSELL *)A->data;
402d1451d4SHong Zhang 
412d1451d4SHong Zhang   PetscFunctionBegin;
422d1451d4SHong Zhang   if (A->offloadmask == PETSC_OFFLOAD_UNALLOCATED || A->offloadmask == PETSC_OFFLOAD_CPU) {
432d1451d4SHong Zhang     PetscCall(PetscLogEventBegin(MAT_CUDACopyToGPU, A, 0, 0, 0));
442d1451d4SHong Zhang     if (A->assembled && A->nonzerostate == cudastruct->nonzerostate) {
452d1451d4SHong Zhang       /* copy values only */
462d1451d4SHong Zhang       PetscCallCUDA(cudaMemcpy(cudastruct->val, a->val, a->sliidx[a->totalslices] * sizeof(MatScalar), cudaMemcpyHostToDevice));
472d1451d4SHong Zhang       PetscCall(PetscLogCpuToGpu(a->sliidx[a->totalslices] * (sizeof(MatScalar))));
482d1451d4SHong Zhang     } else {
492d1451d4SHong Zhang       if (cudastruct->colidx) { PetscCallCUDA(cudaFree(cudastruct->colidx)); }
502d1451d4SHong Zhang       if (cudastruct->val) { PetscCallCUDA(cudaFree(cudastruct->val)); }
512d1451d4SHong Zhang       if (cudastruct->sliidx) { PetscCallCUDA(cudaFree(cudastruct->sliidx)); }
5290d2215bSHong Zhang       if (cudastruct->chunk_slice_map) { PetscCallCUDA(cudaFree(cudastruct->chunk_slice_map)); }
5390d2215bSHong Zhang       cudastruct->maxallocmat  = a->maxallocmat;
5490d2215bSHong Zhang       cudastruct->totalentries = a->sliidx[a->totalslices];
5590d2215bSHong Zhang       cudastruct->totalslices  = a->totalslices;
5690d2215bSHong Zhang       cudastruct->totalchunks  = a->totalchunks;
57*8df136f9SHong Zhang       PetscCallCUDA(cudaMalloc((void **)&(cudastruct->colidx), a->maxallocmat * sizeof(*(cudastruct->colidx))));
58*8df136f9SHong Zhang       PetscCallCUDA(cudaMalloc((void **)&(cudastruct->val), a->maxallocmat * sizeof(*(cudastruct->val))));
592d1451d4SHong Zhang       /* copy values, nz or maxallocmat? */
60*8df136f9SHong Zhang       PetscCallCUDA(cudaMemcpy(cudastruct->colidx, a->colidx, a->sliidx[a->totalslices] * sizeof(*(a->colidx)), cudaMemcpyHostToDevice));
61*8df136f9SHong Zhang       PetscCallCUDA(cudaMemcpy(cudastruct->val, a->val, a->sliidx[a->totalslices] * sizeof(*(a->val)), cudaMemcpyHostToDevice));
6207e43b41SHong Zhang 
63*8df136f9SHong Zhang       PetscCallCUDA(cudaMalloc((void **)&(cudastruct->sliidx), (a->totalslices + 1) * sizeof(*(cudastruct->sliidx))));
64*8df136f9SHong Zhang       PetscCallCUDA(cudaMemcpy(cudastruct->sliidx, a->sliidx, (a->totalslices + 1) * sizeof(*(a->sliidx)), cudaMemcpyHostToDevice));
65*8df136f9SHong Zhang       PetscCallCUDA(cudaMalloc((void **)&(cudastruct->chunk_slice_map), a->totalchunks * sizeof(*(cudastruct->chunk_slice_map))));
66*8df136f9SHong Zhang       PetscCallCUDA(cudaMemcpy(cudastruct->chunk_slice_map, a->chunk_slice_map, a->totalchunks * sizeof(*(a->chunk_slice_map)), cudaMemcpyHostToDevice));
6790d2215bSHong Zhang       PetscCall(PetscLogCpuToGpu(a->sliidx[a->totalslices] * (sizeof(MatScalar) + sizeof(PetscInt)) + (a->totalslices + 1 + a->totalchunks) * sizeof(PetscInt)));
682d1451d4SHong Zhang     }
692d1451d4SHong Zhang     PetscCallCUDA(WaitForCUDA());
702d1451d4SHong Zhang     PetscCall(PetscLogEventEnd(MAT_CUDACopyToGPU, A, 0, 0, 0));
712d1451d4SHong Zhang     A->offloadmask = PETSC_OFFLOAD_BOTH;
722d1451d4SHong Zhang   }
732d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
742d1451d4SHong Zhang }
752d1451d4SHong Zhang 
764e58db63SHong Zhang __global__ void matmult_seqsell_basic_kernel(PetscInt nrows, PetscInt sliceheight, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y)
772d1451d4SHong Zhang {
782d1451d4SHong Zhang   PetscInt  i, row, slice_id, row_in_slice;
792d1451d4SHong Zhang   MatScalar sum;
802d1451d4SHong Zhang   /* one thread per row. */
812d1451d4SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
822d1451d4SHong Zhang   if (row < nrows) {
834e58db63SHong Zhang     slice_id     = row / sliceheight;
844e58db63SHong Zhang     row_in_slice = row % sliceheight;
852d1451d4SHong Zhang     sum          = 0.0;
864e58db63SHong Zhang     for (i = sliidx[slice_id] + row_in_slice; i < sliidx[slice_id + 1]; i += sliceheight) sum += aval[i] * x[acolidx[i]];
872d1451d4SHong Zhang     y[row] = sum;
882d1451d4SHong Zhang   }
892d1451d4SHong Zhang }
902d1451d4SHong Zhang 
914e58db63SHong Zhang __global__ void matmultadd_seqsell_basic_kernel(PetscInt nrows, PetscInt sliceheight, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z)
922d1451d4SHong Zhang {
932d1451d4SHong Zhang   PetscInt  i, row, slice_id, row_in_slice;
942d1451d4SHong Zhang   MatScalar sum;
952d1451d4SHong Zhang   /* one thread per row. */
962d1451d4SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
972d1451d4SHong Zhang   if (row < nrows) {
984e58db63SHong Zhang     slice_id     = row / sliceheight;
994e58db63SHong Zhang     row_in_slice = row % sliceheight;
1002d1451d4SHong Zhang     sum          = 0.0;
1014e58db63SHong Zhang     for (i = sliidx[slice_id] + row_in_slice; i < sliidx[slice_id + 1]; i += sliceheight) sum += aval[i] * x[acolidx[i]];
1022d1451d4SHong Zhang     z[row] = y[row] + sum;
1032d1451d4SHong Zhang   }
1042d1451d4SHong Zhang }
10507e43b41SHong Zhang 
1068711c661SHong Zhang #if !defined(PETSC_USE_COMPLEX)
10707e43b41SHong Zhang /* use 1 block per slice, suitable for large slice width */
10807e43b41SHong Zhang template <int BLOCKY>
1094e58db63SHong Zhang __global__ void matmult_seqsell_tiled_kernel9(PetscInt nrows, PetscInt sliceheight, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y)
11007e43b41SHong Zhang {
1114e58db63SHong Zhang   __shared__ MatScalar shared[32][BLOCKY];
1124e58db63SHong Zhang   PetscInt             i, row, slice_id = blockIdx.x;
1134e58db63SHong Zhang   int                  tid = threadIdx.x + threadIdx.y * 32;
11407e43b41SHong Zhang   /* transposed index */
11507e43b41SHong Zhang   int         tidx = tid % BLOCKY;
11607e43b41SHong Zhang   int         tidy = tid / BLOCKY;
11707e43b41SHong Zhang   PetscScalar t    = 0.0;
1184e58db63SHong Zhang 
1194e58db63SHong Zhang   row = slice_id * sliceheight + threadIdx.x % sliceheight;
12007e43b41SHong Zhang   if (row < nrows) {
1214e58db63SHong Zhang     for (i = sliidx[slice_id] + threadIdx.x + 32 * threadIdx.y; i < sliidx[slice_id + 1]; i += 32 * BLOCKY) t += aval[i] * x[acolidx[i]];
1222d1451d4SHong Zhang   }
1234e58db63SHong Zhang   #pragma unroll
1244e58db63SHong Zhang   for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); }
12507e43b41SHong Zhang   /* transpose layout to reduce each row using warp shfl */
1261f0d1278SHong Zhang   if (threadIdx.x < sliceheight) shared[threadIdx.x][threadIdx.y] = t;
12707e43b41SHong Zhang   __syncthreads();
1281f0d1278SHong Zhang   if (tidy < sliceheight) t = shared[tidy][tidx];
12907e43b41SHong Zhang   #pragma unroll
13007e43b41SHong Zhang   for (int offset = BLOCKY / 2; offset > 0; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset, BLOCKY); }
1314e58db63SHong Zhang   if (tidx == 0 && tidy < sliceheight) { shared[0][tidy] = t; }
13207e43b41SHong Zhang   __syncthreads();
1334e58db63SHong Zhang   if (row < nrows && threadIdx.y == 0 && threadIdx.x < sliceheight) { y[row] = shared[0][threadIdx.x]; }
1342d1451d4SHong Zhang }
1352d1451d4SHong Zhang 
136cca9ff8bSHong Zhang /* use 1 block per slice, suitable for large slice width */
137cca9ff8bSHong Zhang template <int BLOCKY>
138cca9ff8bSHong Zhang __global__ void matmultadd_seqsell_tiled_kernel9(PetscInt nrows, PetscInt sliceheight, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z)
139cca9ff8bSHong Zhang {
140cca9ff8bSHong Zhang   __shared__ MatScalar shared[32][BLOCKY];
141cca9ff8bSHong Zhang   PetscInt             i, row, slice_id = blockIdx.x;
142cca9ff8bSHong Zhang   int                  tid = threadIdx.x + threadIdx.y * 32;
143cca9ff8bSHong Zhang   /* transposed index */
144cca9ff8bSHong Zhang   int         tidx = tid % BLOCKY;
145cca9ff8bSHong Zhang   int         tidy = tid / BLOCKY;
146cca9ff8bSHong Zhang   PetscScalar t    = 0.0;
147cca9ff8bSHong Zhang 
148cca9ff8bSHong Zhang   row = slice_id * sliceheight + threadIdx.x % sliceheight;
149cca9ff8bSHong Zhang   if (row < nrows) {
150cca9ff8bSHong Zhang     for (i = sliidx[slice_id] + threadIdx.x + 32 * threadIdx.y; i < sliidx[slice_id + 1]; i += 32 * BLOCKY) t += aval[i] * x[acolidx[i]];
151cca9ff8bSHong Zhang   }
152cca9ff8bSHong Zhang   #pragma unroll
153cca9ff8bSHong Zhang   for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); }
154cca9ff8bSHong Zhang   /* transpose layout to reduce each row using warp shfl */
1551f0d1278SHong Zhang   if (threadIdx.x < sliceheight) shared[threadIdx.x][threadIdx.y] = t;
156cca9ff8bSHong Zhang   __syncthreads();
1571f0d1278SHong Zhang   if (tidy < sliceheight) t = shared[tidy][tidx];
158cca9ff8bSHong Zhang   #pragma unroll
159cca9ff8bSHong Zhang   for (int offset = BLOCKY / 2; offset > 0; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset, BLOCKY); }
160cca9ff8bSHong Zhang   if (tidx == 0 && tidy < sliceheight) { shared[0][tidy] = t; }
161cca9ff8bSHong Zhang   __syncthreads();
162cca9ff8bSHong Zhang   if (row < nrows && threadIdx.y == 0 && threadIdx.x < sliceheight) { z[row] = y[row] + shared[0][threadIdx.x]; }
163cca9ff8bSHong Zhang }
164cca9ff8bSHong Zhang 
16590d2215bSHong Zhang template <int BLOCKY>
16690d2215bSHong Zhang __device__ __forceinline__ bool segment_scan(PetscInt flag[], MatScalar shared[], PetscScalar *val)
16790d2215bSHong Zhang {
16890d2215bSHong Zhang   bool head = true;
16990d2215bSHong Zhang   #pragma unroll
17090d2215bSHong Zhang   for (int i = 1; i < BLOCKY * 2; i <<= 1) {
17190d2215bSHong Zhang     int halfwarpid                         = threadIdx.y * 2 + threadIdx.x / 16;
17290d2215bSHong Zhang     shared[threadIdx.x + threadIdx.y * 32] = 0;
17390d2215bSHong Zhang     if (halfwarpid >= i && flag[halfwarpid - i] == flag[halfwarpid]) {
17490d2215bSHong Zhang       shared[threadIdx.x + threadIdx.y * 32] = *val;
17590d2215bSHong Zhang       if (i == 1) head = false;
17690d2215bSHong Zhang     }
17790d2215bSHong Zhang     __syncthreads();
17890d2215bSHong Zhang     if (halfwarpid < BLOCKY * 2 - i) *val += shared[threadIdx.x + threadIdx.y * 32 + i * 16];
17990d2215bSHong Zhang     __syncthreads();
18090d2215bSHong Zhang   }
18190d2215bSHong Zhang   return head;
18290d2215bSHong Zhang }
18390d2215bSHong Zhang 
18490d2215bSHong Zhang /* load-balancing version. Chunksize is equal to the number of threads per block */
18590d2215bSHong Zhang template <int BLOCKY>
18690d2215bSHong Zhang __global__ void matmult_seqsell_tiled_kernel8(PetscInt nrows, PetscInt sliceheight, PetscInt chunksperblock, PetscInt totalchunks, const PetscInt *chunk_slice_map, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y)
18790d2215bSHong Zhang {
18890d2215bSHong Zhang   __shared__ MatScalar shared[BLOCKY * 32];
18990d2215bSHong Zhang   PetscInt             gid, row, start_slice, cid;
19090d2215bSHong Zhang   PetscScalar          t = 0.0;
19190d2215bSHong Zhang   /* zero out y */
19290d2215bSHong Zhang   for (int iter = 0; iter < 1 + (nrows - 1) / (gridDim.x * 32 * BLOCKY); iter++) {
19390d2215bSHong Zhang     gid = gridDim.x * 32 * BLOCKY * iter + blockIdx.x * BLOCKY * 32 + threadIdx.y * 32 + threadIdx.x;
19490d2215bSHong Zhang     if (gid < nrows) y[gid] = 0.0;
19590d2215bSHong Zhang   }
19690d2215bSHong Zhang   for (int iter = 0; iter < chunksperblock; iter++) {
19790d2215bSHong Zhang     cid = blockIdx.x * chunksperblock + iter; /* chunk id */
19890d2215bSHong Zhang     if (cid < totalchunks) {
19990d2215bSHong Zhang       start_slice = chunk_slice_map[cid]; /* starting slice at each iteration */
20090d2215bSHong Zhang       gid         = cid * BLOCKY * 32 + threadIdx.y * 32 + threadIdx.x;
20190d2215bSHong Zhang       if ((cid + 1) * BLOCKY * 32 > sliidx[start_slice + 1]) { /* this iteration covers more than one slice */
20290d2215bSHong Zhang         __shared__ PetscInt flag[BLOCKY * 2];
20390d2215bSHong Zhang         bool                write;
20490d2215bSHong Zhang         PetscInt            slice_id = start_slice, totalslices = PetscCeilInt(nrows, sliceheight), totalentries = sliidx[totalslices];
20590d2215bSHong Zhang         /* find out the slice that this element belongs to */
20690d2215bSHong Zhang         while (gid < totalentries && gid >= sliidx[slice_id + 1]) slice_id++;
20790d2215bSHong Zhang         if (threadIdx.x % 16 == 0) flag[threadIdx.y * 2 + threadIdx.x / 16] = slice_id;
20890d2215bSHong Zhang         row = slice_id * sliceheight + threadIdx.x % sliceheight;
20990d2215bSHong Zhang         if (row < nrows && gid < totalentries) t = aval[gid] * x[acolidx[gid]];
21090d2215bSHong Zhang         __syncthreads();
21190d2215bSHong Zhang         write = segment_scan<BLOCKY>(flag, shared, &t);
21290d2215bSHong Zhang         if (row < nrows && gid < totalentries && write) atomicAdd(&y[row], t);
21390d2215bSHong Zhang         t = 0.0;
21490d2215bSHong Zhang       } else { /* this iteration covers only one slice */
21590d2215bSHong Zhang         row = start_slice * sliceheight + threadIdx.x % sliceheight;
21690d2215bSHong Zhang         if (row < nrows) t += aval[gid] * x[acolidx[gid]];
21790d2215bSHong Zhang         if (iter == chunksperblock - 1 || (cid + 2) * BLOCKY * 32 > sliidx[start_slice + 1]) { /* last iteration or next iteration covers more than one slice */
21890d2215bSHong Zhang           int tid = threadIdx.x + threadIdx.y * 32, tidx = tid % BLOCKY, tidy = tid / BLOCKY;
21990d2215bSHong Zhang   /* reduction and write to output vector */
22090d2215bSHong Zhang   #pragma unroll
22190d2215bSHong Zhang           for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); }
22290d2215bSHong Zhang           /* transpose layout to reduce each row using warp shfl */
22390d2215bSHong Zhang           if (threadIdx.x < sliceheight) shared[threadIdx.x * BLOCKY + threadIdx.y] = t; /* shared[threadIdx.x][threadIdx.y] = t */
22490d2215bSHong Zhang           __syncthreads();
22590d2215bSHong Zhang           if (tidy < sliceheight) t = shared[tidy * BLOCKY + tidx]; /* shared[tidy][tidx] */
22690d2215bSHong Zhang   #pragma unroll
22790d2215bSHong Zhang           for (int offset = BLOCKY / 2; offset > 0; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset, BLOCKY); }
22890d2215bSHong Zhang           if (tidx == 0 && tidy < sliceheight) { shared[tidy] = t; /* shared[0][tidy] = t */ }
22990d2215bSHong Zhang           __syncthreads();
23090d2215bSHong Zhang           if (row < nrows && threadIdx.y == 0 && threadIdx.x < sliceheight) atomicAdd(&y[row], shared[threadIdx.x]); /* shared[0][threadIdx.x] */
23190d2215bSHong Zhang           t = 0.0;
23290d2215bSHong Zhang         }
23390d2215bSHong Zhang       }
23490d2215bSHong Zhang     }
23590d2215bSHong Zhang   }
23690d2215bSHong Zhang }
23790d2215bSHong Zhang 
23890d2215bSHong Zhang /* load-balancing version. Chunksize is equal to the number of threads per block */
23990d2215bSHong Zhang template <int BLOCKY>
24090d2215bSHong Zhang __global__ void matmultadd_seqsell_tiled_kernel8(PetscInt nrows, PetscInt sliceheight, PetscInt chunksperblock, PetscInt totalchunks, const PetscInt *chunk_slice_map, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z)
24190d2215bSHong Zhang {
24290d2215bSHong Zhang   __shared__ MatScalar shared[BLOCKY * 32];
24390d2215bSHong Zhang   PetscInt             gid, row, start_slice, cid;
24490d2215bSHong Zhang   PetscScalar          t = 0.0;
24590d2215bSHong Zhang   /* copy y to z */
24690d2215bSHong Zhang   for (int iter = 0; iter < 1 + (nrows - 1) / (gridDim.x * 32 * BLOCKY); iter++) {
24790d2215bSHong Zhang     gid = gridDim.x * 32 * BLOCKY * iter + blockIdx.x * BLOCKY * 32 + threadIdx.y * 32 + threadIdx.x;
24890d2215bSHong Zhang     if (gid < nrows) z[gid] = y[gid];
24990d2215bSHong Zhang   }
25090d2215bSHong Zhang   for (int iter = 0; iter < chunksperblock; iter++) {
25190d2215bSHong Zhang     cid = blockIdx.x * chunksperblock + iter; /* chunk id */
25290d2215bSHong Zhang     if (cid < totalchunks) {
25390d2215bSHong Zhang       start_slice = chunk_slice_map[cid]; /* starting slice at each iteration */
25490d2215bSHong Zhang       gid         = cid * BLOCKY * 32 + threadIdx.y * 32 + threadIdx.x;
25590d2215bSHong Zhang       if ((cid + 1) * BLOCKY * 32 > sliidx[start_slice + 1]) { /* this iteration covers more than one slice */
25690d2215bSHong Zhang         __shared__ PetscInt flag[BLOCKY * 2];
25790d2215bSHong Zhang         bool                write;
25890d2215bSHong Zhang         PetscInt            slice_id = start_slice, totalslices = PetscCeilInt(nrows, sliceheight), totalentries = sliidx[totalslices];
25990d2215bSHong Zhang         /* find out the slice that this element belongs to */
26090d2215bSHong Zhang         while (gid < totalentries && gid >= sliidx[slice_id + 1]) slice_id++;
26190d2215bSHong Zhang         if (threadIdx.x % 16 == 0) flag[threadIdx.y * 2 + threadIdx.x / 16] = slice_id;
26290d2215bSHong Zhang         row = slice_id * sliceheight + threadIdx.x % sliceheight;
26390d2215bSHong Zhang         if (row < nrows && gid < totalentries) t = aval[gid] * x[acolidx[gid]];
26490d2215bSHong Zhang         __syncthreads();
26590d2215bSHong Zhang         write = segment_scan<BLOCKY>(flag, shared, &t);
26690d2215bSHong Zhang         if (row < nrows && gid < totalentries && write) atomicAdd(&z[row], t);
26790d2215bSHong Zhang         t = 0.0;
26890d2215bSHong Zhang       } else { /* this iteration covers only one slice */
26990d2215bSHong Zhang         row = start_slice * sliceheight + threadIdx.x % sliceheight;
27090d2215bSHong Zhang         if (row < nrows) t += aval[gid] * x[acolidx[gid]];
27190d2215bSHong Zhang         if (iter == chunksperblock - 1 || (cid + 2) * BLOCKY * 32 > sliidx[start_slice + 1]) { /* last iteration or next iteration covers more than one slice */
27290d2215bSHong Zhang           int tid = threadIdx.x + threadIdx.y * 32, tidx = tid % BLOCKY, tidy = tid / BLOCKY;
27390d2215bSHong Zhang   /* reduction and write to output vector */
27490d2215bSHong Zhang   #pragma unroll
27590d2215bSHong Zhang           for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); }
27690d2215bSHong Zhang           /* transpose layout to reduce each row using warp shfl */
27790d2215bSHong Zhang           if (threadIdx.x < sliceheight) shared[threadIdx.x * BLOCKY + threadIdx.y] = t; /* shared[threadIdx.x][threadIdx.y] = t */
27890d2215bSHong Zhang           __syncthreads();
27990d2215bSHong Zhang           if (tidy < sliceheight) t = shared[tidy * BLOCKY + tidx]; /* shared[tidy][tidx] */
28090d2215bSHong Zhang   #pragma unroll
28190d2215bSHong Zhang           for (int offset = BLOCKY / 2; offset > 0; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset, BLOCKY); }
28290d2215bSHong Zhang           if (tidx == 0 && tidy < sliceheight) { shared[tidy] = t; /* shared[0][tidy] = t */ }
28390d2215bSHong Zhang           __syncthreads();
28490d2215bSHong Zhang           if (row < nrows && threadIdx.y == 0 && threadIdx.x < sliceheight) atomicAdd(&z[row], shared[threadIdx.x]); /* shared[0][threadIdx.x] */
28590d2215bSHong Zhang           t = 0.0;
28690d2215bSHong Zhang         }
28790d2215bSHong Zhang       }
28890d2215bSHong Zhang     }
28990d2215bSHong Zhang   }
29090d2215bSHong Zhang }
29190d2215bSHong Zhang 
29207e43b41SHong Zhang /* use 1 warp per slice, suitable for small slice width */
2934e58db63SHong Zhang __global__ void matmult_seqsell_tiled_kernel7(PetscInt nrows, PetscInt sliceheight, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y)
2942d1451d4SHong Zhang {
29507e43b41SHong Zhang   PetscInt i, row, slice_id;
29607e43b41SHong Zhang   slice_id = blockIdx.x * blockDim.y + threadIdx.y;
2974e58db63SHong Zhang   row      = slice_id * sliceheight + threadIdx.x % sliceheight;
29807e43b41SHong Zhang   double t = 0.0;
29907e43b41SHong Zhang   if (row < nrows) {
30007e43b41SHong Zhang     for (i = sliidx[slice_id] + threadIdx.x; i < sliidx[slice_id + 1]; i += 32) t += aval[i] * x[acolidx[i]];
30107e43b41SHong Zhang   }
3024e58db63SHong Zhang   #pragma unroll
3034e58db63SHong Zhang   for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); }
3044e58db63SHong Zhang   if (row < nrows && threadIdx.x < sliceheight) { y[row] = t; }
30507e43b41SHong Zhang }
30607e43b41SHong Zhang 
307cca9ff8bSHong Zhang /* use 1 warp per slice, suitable for small slice width */
308cca9ff8bSHong Zhang __global__ void matmultadd_seqsell_tiled_kernel7(PetscInt nrows, PetscInt sliceheight, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z)
309cca9ff8bSHong Zhang {
310cca9ff8bSHong Zhang   PetscInt i, row, slice_id;
311cca9ff8bSHong Zhang   slice_id = blockIdx.x * blockDim.y + threadIdx.y;
312cca9ff8bSHong Zhang   row      = slice_id * sliceheight + threadIdx.x % sliceheight;
313cca9ff8bSHong Zhang   double t = 0.0;
314cca9ff8bSHong Zhang   if (row < nrows) {
315cca9ff8bSHong Zhang     for (i = sliidx[slice_id] + threadIdx.x; i < sliidx[slice_id + 1]; i += 32) t += aval[i] * x[acolidx[i]];
316cca9ff8bSHong Zhang   }
317cca9ff8bSHong Zhang   #pragma unroll
318cca9ff8bSHong Zhang   for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); }
319cca9ff8bSHong Zhang   if (row < nrows && threadIdx.x < sliceheight) { z[row] = y[row] + t; }
320cca9ff8bSHong Zhang }
3218711c661SHong Zhang #endif
322cca9ff8bSHong Zhang 
323a9dd396cSHong Zhang /***********  Kernel 2-6  are tied to slice height 16. They are kept only for performance comparison  **********/
324a9dd396cSHong Zhang 
325a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel6(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y)
32607e43b41SHong Zhang {
32707e43b41SHong Zhang   __shared__ MatScalar shared[512];
3282d1451d4SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
32907e43b41SHong Zhang   /* multiple threads per row. */
3302d1451d4SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
3312d1451d4SHong Zhang   if (row < nrows) {
3322d1451d4SHong Zhang     slice_id     = row / SLICE_HEIGHT;
3332d1451d4SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
3342d1451d4SHong Zhang 
3352d1451d4SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
3362d1451d4SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
33707e43b41SHong Zhang     __syncthreads();
33807e43b41SHong Zhang     if (threadIdx.y < 16) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 16) * blockDim.x + threadIdx.x]; }
33907e43b41SHong Zhang     __syncthreads();
34007e43b41SHong Zhang     if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; }
3412d1451d4SHong Zhang     __syncthreads();
3422d1451d4SHong Zhang     if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; }
3432d1451d4SHong Zhang     __syncthreads();
3442d1451d4SHong Zhang     if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; }
3452d1451d4SHong Zhang     __syncthreads();
3462d1451d4SHong Zhang     if (threadIdx.y < 1) {
34707e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
3482d1451d4SHong Zhang       y[row] = shared[threadIdx.x];
3492d1451d4SHong Zhang     }
3502d1451d4SHong Zhang   }
3512d1451d4SHong Zhang }
3522d1451d4SHong Zhang 
353a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel5(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y)
3542d1451d4SHong Zhang {
35507e43b41SHong Zhang   __shared__ MatScalar shared[512];
3562d1451d4SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
35707e43b41SHong Zhang   /* multiple threads per row. */
3582d1451d4SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
3592d1451d4SHong Zhang   if (row < nrows) {
3602d1451d4SHong Zhang     slice_id     = row / SLICE_HEIGHT;
3612d1451d4SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
3622d1451d4SHong Zhang 
3632d1451d4SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
3642d1451d4SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
36507e43b41SHong Zhang     __syncthreads();
36607e43b41SHong Zhang     if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; }
3672d1451d4SHong Zhang     __syncthreads();
3682d1451d4SHong Zhang     if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; }
3692d1451d4SHong Zhang     __syncthreads();
3702d1451d4SHong Zhang     if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; }
3712d1451d4SHong Zhang     __syncthreads();
3722d1451d4SHong Zhang     if (threadIdx.y < 1) {
37307e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
37407e43b41SHong Zhang       y[row] = shared[threadIdx.x];
37507e43b41SHong Zhang     }
37607e43b41SHong Zhang   }
37707e43b41SHong Zhang }
37807e43b41SHong Zhang 
379a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel4(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y)
38007e43b41SHong Zhang {
38107e43b41SHong Zhang   __shared__ MatScalar shared[512];
38207e43b41SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
38307e43b41SHong Zhang   /* multiple threads per row. */
38407e43b41SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
38507e43b41SHong Zhang   if (row < nrows) {
38607e43b41SHong Zhang     slice_id     = row / SLICE_HEIGHT;
38707e43b41SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
38807e43b41SHong Zhang 
38907e43b41SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
39007e43b41SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
39107e43b41SHong Zhang     __syncthreads();
39207e43b41SHong Zhang     if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; }
39307e43b41SHong Zhang     __syncthreads();
39407e43b41SHong Zhang     if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; }
39507e43b41SHong Zhang     __syncthreads();
39607e43b41SHong Zhang     if (threadIdx.y < 1) {
39707e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
39807e43b41SHong Zhang       y[row] = shared[threadIdx.x];
39907e43b41SHong Zhang     }
40007e43b41SHong Zhang   }
40107e43b41SHong Zhang }
40207e43b41SHong Zhang 
403a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel3(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y)
40407e43b41SHong Zhang {
40507e43b41SHong Zhang   __shared__ MatScalar shared[512];
40607e43b41SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
40707e43b41SHong Zhang   /* multiple threads per row. */
40807e43b41SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
40907e43b41SHong Zhang   if (row < nrows) {
41007e43b41SHong Zhang     slice_id     = row / SLICE_HEIGHT;
41107e43b41SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
41207e43b41SHong Zhang 
41307e43b41SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
41407e43b41SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
41507e43b41SHong Zhang     __syncthreads();
41607e43b41SHong Zhang     if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; }
41707e43b41SHong Zhang     __syncthreads();
41807e43b41SHong Zhang     if (threadIdx.y < 1) {
41907e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
42007e43b41SHong Zhang       y[row] = shared[threadIdx.x];
42107e43b41SHong Zhang     }
42207e43b41SHong Zhang   }
42307e43b41SHong Zhang }
42407e43b41SHong Zhang 
425a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel2(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y)
42607e43b41SHong Zhang {
42707e43b41SHong Zhang   __shared__ MatScalar shared[512];
42807e43b41SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
42907e43b41SHong Zhang   /* multiple threads per row. */
43007e43b41SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
43107e43b41SHong Zhang   if (row < nrows) {
43207e43b41SHong Zhang     slice_id     = row / SLICE_HEIGHT;
43307e43b41SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
43407e43b41SHong Zhang 
43507e43b41SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
43607e43b41SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
43707e43b41SHong Zhang     __syncthreads();
43807e43b41SHong Zhang     if (threadIdx.y < 1) {
43907e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
44007e43b41SHong Zhang       y[row] = shared[threadIdx.x];
44107e43b41SHong Zhang     }
44207e43b41SHong Zhang   }
44307e43b41SHong Zhang }
44407e43b41SHong Zhang 
445a9dd396cSHong Zhang __global__ void matmultadd_seqsell_tiled_kernel6(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z)
44607e43b41SHong Zhang {
44707e43b41SHong Zhang   __shared__ MatScalar shared[512];
44807e43b41SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
44907e43b41SHong Zhang   /* multiple threads per row. */
45007e43b41SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
45107e43b41SHong Zhang   if (row < nrows) {
45207e43b41SHong Zhang     slice_id     = row / SLICE_HEIGHT;
45307e43b41SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
45407e43b41SHong Zhang 
45507e43b41SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
45607e43b41SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
45707e43b41SHong Zhang     __syncthreads();
45807e43b41SHong Zhang     if (threadIdx.y < 16) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 16) * blockDim.x + threadIdx.x]; }
45907e43b41SHong Zhang     __syncthreads();
46007e43b41SHong Zhang     if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; }
46107e43b41SHong Zhang     __syncthreads();
46207e43b41SHong Zhang     if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; }
46307e43b41SHong Zhang     __syncthreads();
46407e43b41SHong Zhang     if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; }
46507e43b41SHong Zhang     __syncthreads();
46607e43b41SHong Zhang     if (threadIdx.y < 1) {
46707e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
4682d1451d4SHong Zhang       z[row] = y[row] + shared[threadIdx.x];
4692d1451d4SHong Zhang     }
4702d1451d4SHong Zhang   }
4712d1451d4SHong Zhang }
47207e43b41SHong Zhang 
473a9dd396cSHong Zhang __global__ void matmultadd_seqsell_tiled_kernel5(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z)
47407e43b41SHong Zhang {
47507e43b41SHong Zhang   __shared__ MatScalar shared[512];
47607e43b41SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
47707e43b41SHong Zhang   /* multiple threads per row. */
47807e43b41SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
47907e43b41SHong Zhang   if (row < nrows) {
48007e43b41SHong Zhang     slice_id     = row / SLICE_HEIGHT;
48107e43b41SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
48207e43b41SHong Zhang 
48307e43b41SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
48407e43b41SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
48507e43b41SHong Zhang     __syncthreads();
48607e43b41SHong Zhang     if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; }
48707e43b41SHong Zhang     __syncthreads();
48807e43b41SHong Zhang     if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; }
48907e43b41SHong Zhang     __syncthreads();
49007e43b41SHong Zhang     if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; }
49107e43b41SHong Zhang     __syncthreads();
49207e43b41SHong Zhang     if (threadIdx.y < 1) {
49307e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
49407e43b41SHong Zhang       z[row] = y[row] + shared[threadIdx.x];
49507e43b41SHong Zhang     }
49607e43b41SHong Zhang   }
49707e43b41SHong Zhang }
49807e43b41SHong Zhang 
499a9dd396cSHong Zhang __global__ void matmultadd_seqsell_tiled_kernel4(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z)
50007e43b41SHong Zhang {
50107e43b41SHong Zhang   __shared__ MatScalar shared[512];
50207e43b41SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
50307e43b41SHong Zhang   /* multiple threads per row. */
50407e43b41SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
50507e43b41SHong Zhang   if (row < nrows) {
50607e43b41SHong Zhang     slice_id     = row / SLICE_HEIGHT;
50707e43b41SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
50807e43b41SHong Zhang 
50907e43b41SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
51007e43b41SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
51107e43b41SHong Zhang     __syncthreads();
51207e43b41SHong Zhang     if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; }
51307e43b41SHong Zhang     __syncthreads();
51407e43b41SHong Zhang     if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; }
51507e43b41SHong Zhang     __syncthreads();
51607e43b41SHong Zhang     if (threadIdx.y < 1) {
51707e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
51807e43b41SHong Zhang       z[row] = y[row] + shared[threadIdx.x];
51907e43b41SHong Zhang     }
52007e43b41SHong Zhang   }
52107e43b41SHong Zhang }
52207e43b41SHong Zhang 
523a9dd396cSHong Zhang __global__ void matmultadd_seqsell_tiled_kernel3(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z)
52407e43b41SHong Zhang {
52507e43b41SHong Zhang   __shared__ MatScalar shared[512];
52607e43b41SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
52707e43b41SHong Zhang   /* multiple threads per row. */
52807e43b41SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
52907e43b41SHong Zhang   if (row < nrows) {
53007e43b41SHong Zhang     slice_id     = row / SLICE_HEIGHT;
53107e43b41SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
53207e43b41SHong Zhang 
53307e43b41SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
53407e43b41SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
53507e43b41SHong Zhang     __syncthreads();
53607e43b41SHong Zhang     if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; }
53707e43b41SHong Zhang     __syncthreads();
53807e43b41SHong Zhang     if (threadIdx.y < 1) {
53907e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
54007e43b41SHong Zhang       z[row] = y[row] + shared[threadIdx.x];
54107e43b41SHong Zhang     }
54207e43b41SHong Zhang   }
54307e43b41SHong Zhang }
54407e43b41SHong Zhang 
545a9dd396cSHong Zhang __global__ void matmultadd_seqsell_tiled_kernel2(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z)
54607e43b41SHong Zhang {
54707e43b41SHong Zhang   __shared__ MatScalar shared[512];
54807e43b41SHong Zhang   PetscInt             i, row, slice_id, row_in_slice;
54907e43b41SHong Zhang   /* multiple threads per row. */
55007e43b41SHong Zhang   row = blockIdx.x * blockDim.x + threadIdx.x;
55107e43b41SHong Zhang   if (row < nrows) {
55207e43b41SHong Zhang     slice_id     = row / SLICE_HEIGHT;
55307e43b41SHong Zhang     row_in_slice = row % SLICE_HEIGHT;
55407e43b41SHong Zhang 
55507e43b41SHong Zhang     shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0;
55607e43b41SHong Zhang     for (i = sliidx[slice_id] + row_in_slice + SLICE_HEIGHT * threadIdx.y; i < sliidx[slice_id + 1]; i += SLICE_HEIGHT * blockDim.y) shared[threadIdx.y * blockDim.x + threadIdx.x] += aval[i] * x[acolidx[i]];
55707e43b41SHong Zhang     __syncthreads();
55807e43b41SHong Zhang     if (threadIdx.y < 1) {
55907e43b41SHong Zhang       shared[threadIdx.x] += shared[blockDim.x + threadIdx.x];
56007e43b41SHong Zhang       z[row] = y[row] + shared[threadIdx.x];
56107e43b41SHong Zhang     }
56207e43b41SHong Zhang   }
5632d1451d4SHong Zhang }
5642d1451d4SHong Zhang 
5652d1451d4SHong Zhang PetscErrorCode MatMult_SeqSELLCUDA(Mat A, Vec xx, Vec yy)
5662d1451d4SHong Zhang {
5672d1451d4SHong Zhang   Mat_SeqSELL       *a          = (Mat_SeqSELL *)A->data;
5682d1451d4SHong Zhang   Mat_SeqSELLCUDA   *cudastruct = (Mat_SeqSELLCUDA *)A->spptr;
5692d1451d4SHong Zhang   PetscScalar       *y;
5702d1451d4SHong Zhang   const PetscScalar *x;
571a9dd396cSHong Zhang   PetscInt           nrows = A->rmap->n, sliceheight = a->sliceheight;
57290d2215bSHong Zhang   PetscInt           chunksperblock, nchunks, *chunk_slice_map;
5732d1451d4SHong Zhang   MatScalar         *aval;
5742d1451d4SHong Zhang   PetscInt          *acolidx;
5752d1451d4SHong Zhang   PetscInt          *sliidx;
57607e43b41SHong Zhang   PetscInt           nblocks, blocksize = 512; /* blocksize must be multiple of SLICE_HEIGHT*32 */
57707e43b41SHong Zhang   dim3               block2(256, 2), block4(128, 4), block8(64, 8), block16(32, 16), block32(16, 32);
578b921024eSHong Zhang   PetscReal          maxoveravg;
5792d1451d4SHong Zhang 
5802d1451d4SHong Zhang   PetscFunctionBegin;
5814e58db63SHong Zhang   PetscCheck(32 % sliceheight == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "The kernel requires a slice height be a divisor of 32, but the input matrix has a slice height of %" PetscInt_FMT, sliceheight);
58290d2215bSHong Zhang   PetscCheck(!(cudastruct->kernelchoice >= 2 && cudastruct->kernelchoice <= 6 && sliceheight != SLICE_HEIGHT), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Kernel choices {2-6} requires the slice height of the matrix be 16, but the current slice height is %" PetscInt_FMT, sliceheight);
5832d1451d4SHong Zhang   PetscCall(MatSeqSELLCUDACopyToGPU(A));
5842d1451d4SHong Zhang   /* cudastruct may not be available until MatSeqSELLCUDACopyToGPU() is called */
5852d1451d4SHong Zhang   aval    = cudastruct->val;
5862d1451d4SHong Zhang   acolidx = cudastruct->colidx;
5872d1451d4SHong Zhang   sliidx  = cudastruct->sliidx;
5882d1451d4SHong Zhang 
5892d1451d4SHong Zhang   PetscCall(VecCUDAGetArrayRead(xx, &x));
5902d1451d4SHong Zhang   PetscCall(VecCUDAGetArrayWrite(yy, &y));
5912d1451d4SHong Zhang   PetscCall(PetscLogGpuTimeBegin());
59207e43b41SHong Zhang 
59307e43b41SHong Zhang   switch (cudastruct->kernelchoice) {
5948711c661SHong Zhang #if !defined(PETSC_USE_COMPLEX)
59507e43b41SHong Zhang   case 9:
5964e58db63SHong Zhang     nblocks = 1 + (nrows - 1) / sliceheight;
5974e58db63SHong Zhang     if (cudastruct->blocky == 2) {
5984e58db63SHong Zhang       matmult_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
5994e58db63SHong Zhang     } else if (cudastruct->blocky == 4) {
6004e58db63SHong Zhang       matmult_seqsell_tiled_kernel9<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6014e58db63SHong Zhang     } else if (cudastruct->blocky == 8) {
6024e58db63SHong Zhang       matmult_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6034e58db63SHong Zhang     } else if (cudastruct->blocky == 16) {
6044e58db63SHong Zhang       matmult_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6054e58db63SHong Zhang     } else if (cudastruct->blocky == 32) {
6064e58db63SHong Zhang       matmult_seqsell_tiled_kernel9<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
607cca9ff8bSHong Zhang     } else {
608cca9ff8bSHong Zhang       matmult_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
60907e43b41SHong Zhang     }
61007e43b41SHong Zhang     break;
61107e43b41SHong Zhang   case 7:
6124e58db63SHong Zhang     nblocks = 1 + (nrows - 1) / (2 * sliceheight);
6134e58db63SHong Zhang     if (cudastruct->blocky == 2) {
6144e58db63SHong Zhang       matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6154e58db63SHong Zhang     } else if (cudastruct->blocky == 4) {
6164e58db63SHong Zhang       matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6174e58db63SHong Zhang     } else if (cudastruct->blocky == 8) {
6184e58db63SHong Zhang       matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6194e58db63SHong Zhang     } else if (cudastruct->blocky == 16) {
6204e58db63SHong Zhang       matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6214e58db63SHong Zhang     } else if (cudastruct->blocky == 32) {
6224e58db63SHong Zhang       matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6234e58db63SHong Zhang     } else {
6244e58db63SHong Zhang       matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6254e58db63SHong Zhang     }
62607e43b41SHong Zhang     break;
6278711c661SHong Zhang #endif
62807e43b41SHong Zhang   case 6:
62907e43b41SHong Zhang     nblocks = 1 + (nrows - 1) / (blocksize / 32); /* 1 slice per block if blocksize=512 */
630a9dd396cSHong Zhang     matmult_seqsell_tiled_kernel6<<<nblocks, block32>>>(nrows, acolidx, aval, sliidx, x, y);
63107e43b41SHong Zhang     break;
63207e43b41SHong Zhang   case 5:
63307e43b41SHong Zhang     nblocks = 1 + (nrows - 1) / (blocksize / 16); /* 2 slices per block if blocksize=512*/
634a9dd396cSHong Zhang     matmult_seqsell_tiled_kernel5<<<nblocks, block16>>>(nrows, acolidx, aval, sliidx, x, y);
63507e43b41SHong Zhang     break;
63607e43b41SHong Zhang   case 4:
63707e43b41SHong Zhang     nblocks = 1 + (nrows - 1) / (blocksize / 8); /* 4 slices per block if blocksize=512 */
638a9dd396cSHong Zhang     matmult_seqsell_tiled_kernel4<<<nblocks, block8>>>(nrows, acolidx, aval, sliidx, x, y);
63907e43b41SHong Zhang     break;
64007e43b41SHong Zhang   case 3:
64107e43b41SHong Zhang     nblocks = 1 + (nrows - 1) / (blocksize / 4); /* 8 slices per block if blocksize=512 */
642a9dd396cSHong Zhang     matmult_seqsell_tiled_kernel3<<<nblocks, block4>>>(nrows, acolidx, aval, sliidx, x, y);
64307e43b41SHong Zhang     break;
64407e43b41SHong Zhang   case 2: /* 16 slices per block if blocksize=512 */
64507e43b41SHong Zhang     nblocks = 1 + (nrows - 1) / (blocksize / 2);
646a9dd396cSHong Zhang     matmult_seqsell_tiled_kernel2<<<nblocks, block2>>>(nrows, acolidx, aval, sliidx, x, y);
64707e43b41SHong Zhang     break;
64807e43b41SHong Zhang   case 1: /* 32 slices per block if blocksize=512 */
64907e43b41SHong Zhang     nblocks = 1 + (nrows - 1) / blocksize;
6504e58db63SHong Zhang     matmult_seqsell_basic_kernel<<<nblocks, blocksize>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
65107e43b41SHong Zhang     break;
6528711c661SHong Zhang #if !defined(PETSC_USE_COMPLEX)
65307e43b41SHong Zhang   case 0:
65490d2215bSHong Zhang     maxoveravg = a->maxslicewidth / a->avgslicewidth;
65590d2215bSHong Zhang     if (maxoveravg > 12.0 && maxoveravg / nrows > 0.001) { /* important threshold */
65690d2215bSHong Zhang       /* each block handles approximately one slice */
65790d2215bSHong Zhang       PetscInt blocky = a->chunksize / 32;
65890d2215bSHong Zhang       nchunks         = cudastruct->totalchunks;
65990d2215bSHong Zhang       chunksperblock  = cudastruct->chunksperblock ? cudastruct->chunksperblock : 1 + (cudastruct->totalentries / cudastruct->totalslices - 1) / a->chunksize;
66090d2215bSHong Zhang       nblocks         = 1 + (nchunks - 1) / chunksperblock;
66190d2215bSHong Zhang       chunk_slice_map = cudastruct->chunk_slice_map;
66290d2215bSHong Zhang       if (blocky == 2) {
66390d2215bSHong Zhang         matmult_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y);
66490d2215bSHong Zhang       } else if (blocky == 4) {
66590d2215bSHong Zhang         matmult_seqsell_tiled_kernel8<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y);
66690d2215bSHong Zhang       } else if (blocky == 8) {
66790d2215bSHong Zhang         matmult_seqsell_tiled_kernel8<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y);
66890d2215bSHong Zhang       } else if (blocky == 16) {
66990d2215bSHong Zhang         matmult_seqsell_tiled_kernel8<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y);
67090d2215bSHong Zhang       } else if (blocky == 32) {
67190d2215bSHong Zhang         matmult_seqsell_tiled_kernel8<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y);
67290d2215bSHong Zhang       } else {
67390d2215bSHong Zhang         matmult_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y);
67490d2215bSHong Zhang       }
6752d1451d4SHong Zhang     } else {
676cca9ff8bSHong Zhang       PetscInt avgslicesize = sliceheight * a->avgslicewidth;
67790d2215bSHong Zhang       if (avgslicesize <= 432) {
67890d2215bSHong Zhang         if (sliceheight * a->maxslicewidth < 2048 && nrows > 100000) {
679cca9ff8bSHong Zhang           nblocks = 1 + (nrows - 1) / (2 * sliceheight); /* two slices per block */
680cca9ff8bSHong Zhang           matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
68190d2215bSHong Zhang         } else {
682cca9ff8bSHong Zhang           nblocks = 1 + (nrows - 1) / sliceheight;
683cca9ff8bSHong Zhang           matmult_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
68490d2215bSHong Zhang         }
685cca9ff8bSHong Zhang       } else if (avgslicesize <= 2400) {
686cca9ff8bSHong Zhang         nblocks = 1 + (nrows - 1) / sliceheight;
687cca9ff8bSHong Zhang         matmult_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
688cca9ff8bSHong Zhang       } else {
6894e58db63SHong Zhang         nblocks = 1 + (nrows - 1) / sliceheight;
6904e58db63SHong Zhang         matmult_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y);
6912d1451d4SHong Zhang       }
6922d1451d4SHong Zhang     }
693cca9ff8bSHong Zhang     break;
6948711c661SHong Zhang #endif
6958711c661SHong Zhang   default:
6968711c661SHong Zhang     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unsupported kernel choice %" PetscInt_FMT " for MatMult_SeqSELLCUDA.", cudastruct->kernelchoice);
697cca9ff8bSHong Zhang   }
6982d1451d4SHong Zhang   PetscCall(PetscLogGpuTimeEnd());
6992d1451d4SHong Zhang   PetscCall(VecCUDARestoreArrayRead(xx, &x));
7002d1451d4SHong Zhang   PetscCall(VecCUDARestoreArrayWrite(yy, &y));
7012d1451d4SHong Zhang   PetscCall(PetscLogGpuFlops(2.0 * a->nz - a->nonzerorowcnt));
7022d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
7032d1451d4SHong Zhang }
7042d1451d4SHong Zhang 
7052d1451d4SHong Zhang PetscErrorCode MatMultAdd_SeqSELLCUDA(Mat A, Vec xx, Vec yy, Vec zz)
7062d1451d4SHong Zhang {
7072d1451d4SHong Zhang   Mat_SeqSELL       *a          = (Mat_SeqSELL *)A->data;
7082d1451d4SHong Zhang   Mat_SeqSELLCUDA   *cudastruct = (Mat_SeqSELLCUDA *)A->spptr;
7092d1451d4SHong Zhang   PetscScalar       *z;
7102d1451d4SHong Zhang   const PetscScalar *y, *x;
711a9dd396cSHong Zhang   PetscInt           nrows = A->rmap->n, sliceheight = a->sliceheight;
71290d2215bSHong Zhang   PetscInt           chunksperblock, nchunks, *chunk_slice_map;
7132d1451d4SHong Zhang   MatScalar         *aval    = cudastruct->val;
7142d1451d4SHong Zhang   PetscInt          *acolidx = cudastruct->colidx;
7152d1451d4SHong Zhang   PetscInt          *sliidx  = cudastruct->sliidx;
71690d2215bSHong Zhang   PetscReal          maxoveravg;
7172d1451d4SHong Zhang 
7182d1451d4SHong Zhang   PetscFunctionBegin;
71990d2215bSHong Zhang   PetscCheck(32 % sliceheight == 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "The kernel requires a slice height be a divisor of 32, but the input matrix has a slice height of %" PetscInt_FMT, sliceheight);
72090d2215bSHong Zhang   PetscCheck(!(cudastruct->kernelchoice >= 2 && cudastruct->kernelchoice <= 6 && sliceheight != SLICE_HEIGHT), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Kernel choices {2-6} requires the slice height of the matrix be 16, but the current slice height is %" PetscInt_FMT, sliceheight);
7212d1451d4SHong Zhang   PetscCall(MatSeqSELLCUDACopyToGPU(A));
7222d1451d4SHong Zhang   if (a->nz) {
72390d2215bSHong Zhang     PetscInt blocky = cudastruct->blocky, nblocks, blocksize = 512;
72407e43b41SHong Zhang     dim3     block2(256, 2), block4(128, 4), block8(64, 8), block16(32, 16), block32(16, 32);
7252d1451d4SHong Zhang     PetscCall(VecCUDAGetArrayRead(xx, &x));
7262d1451d4SHong Zhang     PetscCall(VecCUDAGetArrayRead(yy, &y));
7272d1451d4SHong Zhang     PetscCall(VecCUDAGetArrayWrite(zz, &z));
7282d1451d4SHong Zhang     PetscCall(PetscLogGpuTimeBegin());
72907e43b41SHong Zhang 
73007e43b41SHong Zhang     switch (cudastruct->kernelchoice) {
7318711c661SHong Zhang #if !defined(PETSC_USE_COMPLEX)
732cca9ff8bSHong Zhang     case 9:
733cca9ff8bSHong Zhang       nblocks = 1 + (nrows - 1) / sliceheight;
73490d2215bSHong Zhang       if (blocky == 2) {
735cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
73690d2215bSHong Zhang       } else if (blocky == 4) {
737cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel9<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
73890d2215bSHong Zhang       } else if (blocky == 8) {
739cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
74090d2215bSHong Zhang       } else if (blocky == 16) {
741cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
74290d2215bSHong Zhang       } else if (blocky == 32) {
743cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel9<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
744cca9ff8bSHong Zhang       } else {
745cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
746cca9ff8bSHong Zhang       }
747cca9ff8bSHong Zhang       break;
74890d2215bSHong Zhang     case 8:
74990d2215bSHong Zhang       /* each block handles approximately one slice */
75090d2215bSHong Zhang       nchunks         = cudastruct->totalchunks;
75190d2215bSHong Zhang       blocky          = a->chunksize / 32;
75290d2215bSHong Zhang       chunksperblock  = cudastruct->chunksperblock ? cudastruct->chunksperblock : 1 + (cudastruct->totalentries / cudastruct->totalslices - 1) / a->chunksize;
75390d2215bSHong Zhang       nblocks         = 1 + (nchunks - 1) / chunksperblock;
75490d2215bSHong Zhang       chunk_slice_map = cudastruct->chunk_slice_map;
75590d2215bSHong Zhang       if (blocky == 2) {
75690d2215bSHong Zhang         matmultadd_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
75790d2215bSHong Zhang       } else if (blocky == 4) {
75890d2215bSHong Zhang         matmultadd_seqsell_tiled_kernel8<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
75990d2215bSHong Zhang       } else if (blocky == 8) {
76090d2215bSHong Zhang         matmultadd_seqsell_tiled_kernel8<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
76190d2215bSHong Zhang       } else if (blocky == 16) {
76290d2215bSHong Zhang         matmultadd_seqsell_tiled_kernel8<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
76390d2215bSHong Zhang       } else if (blocky == 32) {
76490d2215bSHong Zhang         matmultadd_seqsell_tiled_kernel8<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
76590d2215bSHong Zhang       } else {
76690d2215bSHong Zhang         matmultadd_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
76790d2215bSHong Zhang       }
76890d2215bSHong Zhang       break;
769cca9ff8bSHong Zhang     case 7:
770cca9ff8bSHong Zhang       nblocks = 1 + (nrows - 1) / (2 * sliceheight);
77190d2215bSHong Zhang       if (blocky == 2) {
772cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
77390d2215bSHong Zhang       } else if (blocky == 4) {
774cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
77590d2215bSHong Zhang       } else if (blocky == 8) {
776cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
77790d2215bSHong Zhang       } else if (blocky == 16) {
778cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
77990d2215bSHong Zhang       } else if (blocky == 32) {
780cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
781cca9ff8bSHong Zhang       } else {
782cca9ff8bSHong Zhang         matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
783cca9ff8bSHong Zhang       }
784cca9ff8bSHong Zhang       break;
7858711c661SHong Zhang #endif
78607e43b41SHong Zhang     case 6:
78707e43b41SHong Zhang       nblocks = 1 + (nrows - 1) / (blocksize / 32);
788a9dd396cSHong Zhang       matmultadd_seqsell_tiled_kernel6<<<nblocks, block32>>>(nrows, acolidx, aval, sliidx, x, y, z);
78907e43b41SHong Zhang       break;
79007e43b41SHong Zhang     case 5:
79107e43b41SHong Zhang       nblocks = 1 + (nrows - 1) / (blocksize / 16);
792a9dd396cSHong Zhang       matmultadd_seqsell_tiled_kernel5<<<nblocks, block16>>>(nrows, acolidx, aval, sliidx, x, y, z);
79307e43b41SHong Zhang       break;
79407e43b41SHong Zhang     case 4:
79507e43b41SHong Zhang       nblocks = 1 + (nrows - 1) / (blocksize / 8);
796a9dd396cSHong Zhang       matmultadd_seqsell_tiled_kernel4<<<nblocks, block8>>>(nrows, acolidx, aval, sliidx, x, y, z);
79707e43b41SHong Zhang       break;
79807e43b41SHong Zhang     case 3:
79907e43b41SHong Zhang       nblocks = 1 + (nrows - 1) / (blocksize / 4);
800a9dd396cSHong Zhang       matmultadd_seqsell_tiled_kernel3<<<nblocks, block4>>>(nrows, acolidx, aval, sliidx, x, y, z);
80107e43b41SHong Zhang       break;
80207e43b41SHong Zhang     case 2:
80307e43b41SHong Zhang       nblocks = 1 + (nrows - 1) / (blocksize / 2);
804a9dd396cSHong Zhang       matmultadd_seqsell_tiled_kernel2<<<nblocks, block2>>>(nrows, acolidx, aval, sliidx, x, y, z);
80507e43b41SHong Zhang       break;
80607e43b41SHong Zhang     case 1:
80707e43b41SHong Zhang       nblocks = 1 + (nrows - 1) / blocksize;
8084e58db63SHong Zhang       matmultadd_seqsell_basic_kernel<<<nblocks, blocksize>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
80907e43b41SHong Zhang       break;
8108711c661SHong Zhang #if !defined(PETSC_USE_COMPLEX)
811cca9ff8bSHong Zhang     case 0:
81290d2215bSHong Zhang       maxoveravg = a->maxslicewidth / a->avgslicewidth;
81390d2215bSHong Zhang       if (maxoveravg > 12.0 && maxoveravg / nrows > 0.001) { /* important threshold */
81490d2215bSHong Zhang         /* each block handles approximately one slice */
81590d2215bSHong Zhang         nchunks         = cudastruct->totalchunks;
81690d2215bSHong Zhang         blocky          = a->chunksize / 32;
81790d2215bSHong Zhang         chunksperblock  = cudastruct->chunksperblock ? cudastruct->chunksperblock : 1 + (cudastruct->totalentries / cudastruct->totalslices - 1) / a->chunksize;
81890d2215bSHong Zhang         nblocks         = 1 + (nchunks - 1) / chunksperblock;
81990d2215bSHong Zhang         chunk_slice_map = cudastruct->chunk_slice_map;
82090d2215bSHong Zhang         if (blocky == 2) {
82190d2215bSHong Zhang           matmultadd_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
82290d2215bSHong Zhang         } else if (blocky == 4) {
82390d2215bSHong Zhang           matmultadd_seqsell_tiled_kernel8<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
82490d2215bSHong Zhang         } else if (blocky == 8) {
82590d2215bSHong Zhang           matmultadd_seqsell_tiled_kernel8<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
82690d2215bSHong Zhang         } else if (blocky == 16) {
82790d2215bSHong Zhang           matmultadd_seqsell_tiled_kernel8<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
82890d2215bSHong Zhang         } else if (blocky == 32) {
82990d2215bSHong Zhang           matmultadd_seqsell_tiled_kernel8<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
83090d2215bSHong Zhang         } else {
83190d2215bSHong Zhang           matmultadd_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z);
83290d2215bSHong Zhang         }
833cca9ff8bSHong Zhang       } else {
834cca9ff8bSHong Zhang         PetscInt avgslicesize = sliceheight * a->avgslicewidth;
83590d2215bSHong Zhang         if (avgslicesize <= 432) {
83690d2215bSHong Zhang           if (sliceheight * a->maxslicewidth < 2048 && nrows > 100000) {
837cca9ff8bSHong Zhang             nblocks = 1 + (nrows - 1) / (2 * sliceheight); /* two slices per block */
838cca9ff8bSHong Zhang             matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
83990d2215bSHong Zhang           } else {
840cca9ff8bSHong Zhang             nblocks = 1 + (nrows - 1) / sliceheight;
841cca9ff8bSHong Zhang             matmultadd_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
84290d2215bSHong Zhang           }
843cca9ff8bSHong Zhang         } else if (avgslicesize <= 2400) {
844cca9ff8bSHong Zhang           nblocks = 1 + (nrows - 1) / sliceheight;
845cca9ff8bSHong Zhang           matmultadd_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
846cca9ff8bSHong Zhang         } else {
847cca9ff8bSHong Zhang           nblocks = 1 + (nrows - 1) / sliceheight;
848cca9ff8bSHong Zhang           matmultadd_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z);
849cca9ff8bSHong Zhang         }
850cca9ff8bSHong Zhang       }
85107e43b41SHong Zhang       break;
8528711c661SHong Zhang #endif
8538711c661SHong Zhang     default:
8548711c661SHong Zhang       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "unsupported kernel choice %" PetscInt_FMT " for MatMult_SeqSELLCUDA.", cudastruct->kernelchoice);
8552d1451d4SHong Zhang     }
8562d1451d4SHong Zhang     PetscCall(PetscLogGpuTimeEnd());
8572d1451d4SHong Zhang     PetscCall(VecCUDARestoreArrayRead(xx, &x));
8582d1451d4SHong Zhang     PetscCall(VecCUDARestoreArrayRead(yy, &y));
8592d1451d4SHong Zhang     PetscCall(VecCUDARestoreArrayWrite(zz, &z));
8602d1451d4SHong Zhang     PetscCall(PetscLogGpuFlops(2.0 * a->nz));
8612d1451d4SHong Zhang   } else {
8622d1451d4SHong Zhang     PetscCall(VecCopy(yy, zz));
8632d1451d4SHong Zhang   }
8642d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
8652d1451d4SHong Zhang }
8662d1451d4SHong Zhang 
8672d1451d4SHong Zhang static PetscErrorCode MatSetFromOptions_SeqSELLCUDA(Mat A, PetscOptionItems *PetscOptionsObject)
8682d1451d4SHong Zhang {
86907e43b41SHong Zhang   Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr;
8704e58db63SHong Zhang   PetscInt         kernel, blocky;
87107e43b41SHong Zhang   PetscBool        flg;
87207e43b41SHong Zhang 
8732d1451d4SHong Zhang   PetscFunctionBegin;
8742d1451d4SHong Zhang   PetscOptionsHeadBegin(PetscOptionsObject, "SeqSELLCUDA options");
87590d2215bSHong Zhang   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_sell_spmv_cuda_blocky", &blocky, &flg));
87690d2215bSHong Zhang   if (flg) {
87790d2215bSHong Zhang     PetscCheck(blocky == 2 || blocky == 4 || blocky == 8 || blocky == 16 || blocky == 32, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported blocky: %" PetscInt_FMT " it should be in {2,4,8,16,32}", blocky);
87890d2215bSHong Zhang     cudastruct->blocky = blocky;
87990d2215bSHong Zhang   }
88007e43b41SHong Zhang   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_sell_spmv_cuda_kernel", &kernel, &flg));
88107e43b41SHong Zhang   if (flg) {
88207e43b41SHong Zhang     PetscCheck(kernel >= 0 && kernel <= 9, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Wrong kernel choice: %" PetscInt_FMT " it should be in [0,9]", kernel);
88307e43b41SHong Zhang     cudastruct->kernelchoice = kernel;
88490d2215bSHong Zhang     if (kernel == 8) { PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_sell_spmv_cuda_chunksperblock", &cudastruct->chunksperblock, &flg)); }
8854e58db63SHong Zhang   }
8862d1451d4SHong Zhang   PetscOptionsHeadEnd();
8872d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
8882d1451d4SHong Zhang }
8892d1451d4SHong Zhang 
89007e43b41SHong Zhang PETSC_INTERN PetscErrorCode MatAssemblyEnd_SpMV_Preprocessing_Private(Mat A)
89107e43b41SHong Zhang {
89207e43b41SHong Zhang   Mat_SeqSELL *a = (Mat_SeqSELL *)A->data;
89307e43b41SHong Zhang 
89490d2215bSHong Zhang   PetscFunctionBegin;
89507e43b41SHong Zhang   PetscCall(MatSeqSELLGetAvgSliceWidth(A, &a->avgslicewidth));
89607e43b41SHong Zhang   PetscCall(MatSeqSELLGetMaxSliceWidth(A, &a->maxslicewidth));
89707e43b41SHong Zhang   PetscCall(MatSeqSELLGetFillRatio(A, &a->fillratio));
89807e43b41SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
89907e43b41SHong Zhang }
90007e43b41SHong Zhang 
9012d1451d4SHong Zhang static PetscErrorCode MatAssemblyEnd_SeqSELLCUDA(Mat A, MatAssemblyType mode)
9022d1451d4SHong Zhang {
9032d1451d4SHong Zhang   PetscFunctionBegin;
9042d1451d4SHong Zhang   PetscCall(MatAssemblyEnd_SeqSELL(A, mode));
90507e43b41SHong Zhang   PetscCall(MatAssemblyEnd_SpMV_Preprocessing_Private(A));
9062d1451d4SHong Zhang   if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS);
9072d1451d4SHong Zhang   if (A->factortype == MAT_FACTOR_NONE) { PetscCall(MatSeqSELLCUDACopyToGPU(A)); }
9082d1451d4SHong Zhang   A->ops->mult    = MatMult_SeqSELLCUDA;
9092d1451d4SHong Zhang   A->ops->multadd = MatMultAdd_SeqSELLCUDA;
9102d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
9112d1451d4SHong Zhang }
9122d1451d4SHong Zhang 
913*8df136f9SHong Zhang static PetscErrorCode MatZeroEntries_SeqSELLCUDA(Mat A)
914*8df136f9SHong Zhang {
915*8df136f9SHong Zhang   PetscBool    both = PETSC_FALSE;
916*8df136f9SHong Zhang   Mat_SeqSELL *a    = (Mat_SeqSELL *)A->data;
917*8df136f9SHong Zhang 
918*8df136f9SHong Zhang   PetscFunctionBegin;
919*8df136f9SHong Zhang   if (A->factortype == MAT_FACTOR_NONE) {
920*8df136f9SHong Zhang     Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr;
921*8df136f9SHong Zhang     if (cudastruct->val) {
922*8df136f9SHong Zhang       both = PETSC_TRUE;
923*8df136f9SHong Zhang       PetscCallCUDA(cudaMemset(cudastruct->val, 0, a->sliidx[a->totalslices] * sizeof(*(cudastruct->val))));
924*8df136f9SHong Zhang     }
925*8df136f9SHong Zhang   }
926*8df136f9SHong Zhang   PetscCall(PetscArrayzero(a->val, a->sliidx[a->totalslices]));
927*8df136f9SHong Zhang   PetscCall(MatSeqSELLInvalidateDiagonal(A));
928*8df136f9SHong Zhang   if (both) A->offloadmask = PETSC_OFFLOAD_BOTH;
929*8df136f9SHong Zhang   else A->offloadmask = PETSC_OFFLOAD_CPU;
930*8df136f9SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
931*8df136f9SHong Zhang }
932*8df136f9SHong Zhang 
9332d1451d4SHong Zhang static PetscErrorCode MatDestroy_SeqSELLCUDA(Mat A)
9342d1451d4SHong Zhang {
9352d1451d4SHong Zhang   PetscFunctionBegin;
936*8df136f9SHong Zhang   if (A->factortype == MAT_FACTOR_NONE && A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) PetscCall(MatSeqSELLCUDA_Destroy((Mat_SeqSELLCUDA **)&A->spptr));
9372d1451d4SHong Zhang   PetscCall(MatDestroy_SeqSELL(A));
9382d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
9392d1451d4SHong Zhang }
9402d1451d4SHong Zhang 
9412d1451d4SHong Zhang PETSC_INTERN PetscErrorCode MatConvert_SeqSELL_SeqSELLCUDA(Mat);
9422d1451d4SHong Zhang static PetscErrorCode       MatDuplicate_SeqSELLCUDA(Mat A, MatDuplicateOption cpvalues, Mat *B)
9432d1451d4SHong Zhang {
9442d1451d4SHong Zhang   PetscFunctionBegin;
9452d1451d4SHong Zhang   PetscCall(MatDuplicate_SeqSELL(A, cpvalues, B));
9462d1451d4SHong Zhang   PetscCall(MatConvert_SeqSELL_SeqSELLCUDA(*B));
9472d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
9482d1451d4SHong Zhang }
9492d1451d4SHong Zhang 
950*8df136f9SHong Zhang PETSC_INTERN PetscErrorCode MatConvert_SeqSELL_SeqSELLCUDA(Mat B)
9512d1451d4SHong Zhang {
9522d1451d4SHong Zhang   Mat_SeqSELLCUDA *cudastruct;
9532d1451d4SHong Zhang 
9542d1451d4SHong Zhang   PetscFunctionBegin;
9552d1451d4SHong Zhang   PetscCall(PetscFree(B->defaultvectype));
9562d1451d4SHong Zhang   PetscCall(PetscStrallocpy(VECCUDA, &B->defaultvectype));
9572d1451d4SHong Zhang 
9582d1451d4SHong Zhang   if (!B->spptr) {
9592d1451d4SHong Zhang     if (B->factortype == MAT_FACTOR_NONE) {
9602d1451d4SHong Zhang       PetscCall(PetscNew(&cudastruct));
9612d1451d4SHong Zhang       B->spptr = cudastruct;
9622d1451d4SHong Zhang     }
9632d1451d4SHong Zhang   }
9642d1451d4SHong Zhang 
9652d1451d4SHong Zhang   B->ops->assemblyend    = MatAssemblyEnd_SeqSELLCUDA;
9662d1451d4SHong Zhang   B->ops->destroy        = MatDestroy_SeqSELLCUDA;
9672d1451d4SHong Zhang   B->ops->setfromoptions = MatSetFromOptions_SeqSELLCUDA;
9682d1451d4SHong Zhang   B->ops->mult           = MatMult_SeqSELLCUDA;
9692d1451d4SHong Zhang   B->ops->multadd        = MatMultAdd_SeqSELLCUDA;
9702d1451d4SHong Zhang   B->ops->duplicate      = MatDuplicate_SeqSELLCUDA;
971*8df136f9SHong Zhang   B->ops->zeroentries    = MatZeroEntries_SeqSELLCUDA;
9722d1451d4SHong Zhang 
97307e43b41SHong Zhang   /* No need to assemble SeqSELL, but need to do the preprocessing for SpMV */
97407e43b41SHong Zhang   PetscCall(MatAssemblyEnd_SpMV_Preprocessing_Private(B));
97507e43b41SHong Zhang 
9762d1451d4SHong Zhang   PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQSELLCUDA));
9772d1451d4SHong Zhang   B->offloadmask = PETSC_OFFLOAD_UNALLOCATED;
9782d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
9792d1451d4SHong Zhang }
9802d1451d4SHong Zhang 
9812d1451d4SHong Zhang PETSC_EXTERN PetscErrorCode MatCreate_SeqSELLCUDA(Mat B)
9822d1451d4SHong Zhang {
9832d1451d4SHong Zhang   PetscFunctionBegin;
9842d1451d4SHong Zhang   PetscCall(MatCreate_SeqSELL(B));
9852d1451d4SHong Zhang   PetscCall(MatConvert_SeqSELL_SeqSELLCUDA(B));
986*8df136f9SHong Zhang   PetscCall(MatSetFromOptions(B));
9872d1451d4SHong Zhang   PetscFunctionReturn(PETSC_SUCCESS);
9882d1451d4SHong Zhang }
989