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; 572d1451d4SHong Zhang PetscCallCUDA(cudaMalloc((void **)&(cudastruct->colidx), a->maxallocmat * sizeof(PetscInt))); 582d1451d4SHong Zhang PetscCallCUDA(cudaMalloc((void **)&(cudastruct->val), a->maxallocmat * sizeof(MatScalar))); 592d1451d4SHong Zhang /* copy values, nz or maxallocmat? */ 602d1451d4SHong Zhang PetscCallCUDA(cudaMemcpy(cudastruct->colidx, a->colidx, a->sliidx[a->totalslices] * sizeof(PetscInt), cudaMemcpyHostToDevice)); 612d1451d4SHong Zhang PetscCallCUDA(cudaMemcpy(cudastruct->val, a->val, a->sliidx[a->totalslices] * sizeof(MatScalar), cudaMemcpyHostToDevice)); 6207e43b41SHong Zhang 6307e43b41SHong Zhang PetscCallCUDA(cudaMalloc((void **)&(cudastruct->sliidx), (a->totalslices + 1) * sizeof(PetscInt))); 6407e43b41SHong Zhang PetscCallCUDA(cudaMemcpy(cudastruct->sliidx, a->sliidx, (a->totalslices + 1) * sizeof(PetscInt), cudaMemcpyHostToDevice)); 6590d2215bSHong Zhang PetscCallCUDA(cudaMalloc((void **)&(cudastruct->chunk_slice_map), a->totalchunks * sizeof(PetscInt))); 6690d2215bSHong Zhang PetscCallCUDA(cudaMemcpy(cudastruct->chunk_slice_map, a->chunk_slice_map, a->totalchunks * sizeof(PetscInt), 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 10607e43b41SHong Zhang /* use 1 block per slice, suitable for large slice width */ 10707e43b41SHong Zhang template <int BLOCKY> 1084e58db63SHong 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) 10907e43b41SHong Zhang { 1104e58db63SHong Zhang __shared__ MatScalar shared[32][BLOCKY]; 1114e58db63SHong Zhang PetscInt i, row, slice_id = blockIdx.x; 1124e58db63SHong Zhang int tid = threadIdx.x + threadIdx.y * 32; 11307e43b41SHong Zhang /* transposed index */ 11407e43b41SHong Zhang int tidx = tid % BLOCKY; 11507e43b41SHong Zhang int tidy = tid / BLOCKY; 11607e43b41SHong Zhang PetscScalar t = 0.0; 1174e58db63SHong Zhang 1184e58db63SHong Zhang row = slice_id * sliceheight + threadIdx.x % sliceheight; 11907e43b41SHong Zhang if (row < nrows) { 1204e58db63SHong 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]]; 1212d1451d4SHong Zhang } 1224e58db63SHong Zhang #pragma unroll 1234e58db63SHong Zhang for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); } 12407e43b41SHong Zhang /* transpose layout to reduce each row using warp shfl */ 1251f0d1278SHong Zhang if (threadIdx.x < sliceheight) shared[threadIdx.x][threadIdx.y] = t; 12607e43b41SHong Zhang __syncthreads(); 1271f0d1278SHong Zhang if (tidy < sliceheight) t = shared[tidy][tidx]; 12807e43b41SHong Zhang #pragma unroll 12907e43b41SHong Zhang for (int offset = BLOCKY / 2; offset > 0; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset, BLOCKY); } 1304e58db63SHong Zhang if (tidx == 0 && tidy < sliceheight) { shared[0][tidy] = t; } 13107e43b41SHong Zhang __syncthreads(); 1324e58db63SHong Zhang if (row < nrows && threadIdx.y == 0 && threadIdx.x < sliceheight) { y[row] = shared[0][threadIdx.x]; } 1332d1451d4SHong Zhang } 1342d1451d4SHong Zhang 135cca9ff8bSHong Zhang /* use 1 block per slice, suitable for large slice width */ 136cca9ff8bSHong Zhang template <int BLOCKY> 137cca9ff8bSHong 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) 138cca9ff8bSHong Zhang { 139cca9ff8bSHong Zhang __shared__ MatScalar shared[32][BLOCKY]; 140cca9ff8bSHong Zhang PetscInt i, row, slice_id = blockIdx.x; 141cca9ff8bSHong Zhang int tid = threadIdx.x + threadIdx.y * 32; 142cca9ff8bSHong Zhang /* transposed index */ 143cca9ff8bSHong Zhang int tidx = tid % BLOCKY; 144cca9ff8bSHong Zhang int tidy = tid / BLOCKY; 145cca9ff8bSHong Zhang PetscScalar t = 0.0; 146cca9ff8bSHong Zhang 147cca9ff8bSHong Zhang row = slice_id * sliceheight + threadIdx.x % sliceheight; 148cca9ff8bSHong Zhang if (row < nrows) { 149cca9ff8bSHong 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]]; 150cca9ff8bSHong Zhang } 151cca9ff8bSHong Zhang #pragma unroll 152cca9ff8bSHong Zhang for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); } 153cca9ff8bSHong Zhang /* transpose layout to reduce each row using warp shfl */ 1541f0d1278SHong Zhang if (threadIdx.x < sliceheight) shared[threadIdx.x][threadIdx.y] = t; 155cca9ff8bSHong Zhang __syncthreads(); 1561f0d1278SHong Zhang if (tidy < sliceheight) t = shared[tidy][tidx]; 157cca9ff8bSHong Zhang #pragma unroll 158cca9ff8bSHong Zhang for (int offset = BLOCKY / 2; offset > 0; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset, BLOCKY); } 159cca9ff8bSHong Zhang if (tidx == 0 && tidy < sliceheight) { shared[0][tidy] = t; } 160cca9ff8bSHong Zhang __syncthreads(); 161cca9ff8bSHong Zhang if (row < nrows && threadIdx.y == 0 && threadIdx.x < sliceheight) { z[row] = y[row] + shared[0][threadIdx.x]; } 162cca9ff8bSHong Zhang } 163cca9ff8bSHong Zhang 16490d2215bSHong Zhang template <int BLOCKY> 16590d2215bSHong Zhang __device__ __forceinline__ bool segment_scan(PetscInt flag[], MatScalar shared[], PetscScalar *val) 16690d2215bSHong Zhang { 16790d2215bSHong Zhang bool head = true; 16890d2215bSHong Zhang #pragma unroll 16990d2215bSHong Zhang for (int i = 1; i < BLOCKY * 2; i <<= 1) { 17090d2215bSHong Zhang int halfwarpid = threadIdx.y * 2 + threadIdx.x / 16; 17190d2215bSHong Zhang shared[threadIdx.x + threadIdx.y * 32] = 0; 17290d2215bSHong Zhang if (halfwarpid >= i && flag[halfwarpid - i] == flag[halfwarpid]) { 17390d2215bSHong Zhang shared[threadIdx.x + threadIdx.y * 32] = *val; 17490d2215bSHong Zhang if (i == 1) head = false; 17590d2215bSHong Zhang } 17690d2215bSHong Zhang __syncthreads(); 17790d2215bSHong Zhang if (halfwarpid < BLOCKY * 2 - i) *val += shared[threadIdx.x + threadIdx.y * 32 + i * 16]; 17890d2215bSHong Zhang __syncthreads(); 17990d2215bSHong Zhang } 18090d2215bSHong Zhang return head; 18190d2215bSHong Zhang } 18290d2215bSHong Zhang 18390d2215bSHong Zhang /* load-balancing version. Chunksize is equal to the number of threads per block */ 18490d2215bSHong Zhang template <int BLOCKY> 18590d2215bSHong 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) 18690d2215bSHong Zhang { 18790d2215bSHong Zhang __shared__ MatScalar shared[BLOCKY * 32]; 18890d2215bSHong Zhang PetscInt gid, row, start_slice, cid; 18990d2215bSHong Zhang PetscScalar t = 0.0; 19090d2215bSHong Zhang /* zero out y */ 19190d2215bSHong Zhang for (int iter = 0; iter < 1 + (nrows - 1) / (gridDim.x * 32 * BLOCKY); iter++) { 19290d2215bSHong Zhang gid = gridDim.x * 32 * BLOCKY * iter + blockIdx.x * BLOCKY * 32 + threadIdx.y * 32 + threadIdx.x; 19390d2215bSHong Zhang if (gid < nrows) y[gid] = 0.0; 19490d2215bSHong Zhang } 19590d2215bSHong Zhang for (int iter = 0; iter < chunksperblock; iter++) { 19690d2215bSHong Zhang cid = blockIdx.x * chunksperblock + iter; /* chunk id */ 19790d2215bSHong Zhang if (cid < totalchunks) { 19890d2215bSHong Zhang start_slice = chunk_slice_map[cid]; /* starting slice at each iteration */ 19990d2215bSHong Zhang gid = cid * BLOCKY * 32 + threadIdx.y * 32 + threadIdx.x; 20090d2215bSHong Zhang if ((cid + 1) * BLOCKY * 32 > sliidx[start_slice + 1]) { /* this iteration covers more than one slice */ 20190d2215bSHong Zhang __shared__ PetscInt flag[BLOCKY * 2]; 20290d2215bSHong Zhang bool write; 20390d2215bSHong Zhang PetscInt slice_id = start_slice, totalslices = PetscCeilInt(nrows, sliceheight), totalentries = sliidx[totalslices]; 20490d2215bSHong Zhang /* find out the slice that this element belongs to */ 20590d2215bSHong Zhang while (gid < totalentries && gid >= sliidx[slice_id + 1]) slice_id++; 20690d2215bSHong Zhang if (threadIdx.x % 16 == 0) flag[threadIdx.y * 2 + threadIdx.x / 16] = slice_id; 20790d2215bSHong Zhang row = slice_id * sliceheight + threadIdx.x % sliceheight; 20890d2215bSHong Zhang if (row < nrows && gid < totalentries) t = aval[gid] * x[acolidx[gid]]; 20990d2215bSHong Zhang __syncthreads(); 21090d2215bSHong Zhang write = segment_scan<BLOCKY>(flag, shared, &t); 21190d2215bSHong Zhang if (row < nrows && gid < totalentries && write) atomicAdd(&y[row], t); 21290d2215bSHong Zhang t = 0.0; 21390d2215bSHong Zhang } else { /* this iteration covers only one slice */ 21490d2215bSHong Zhang row = start_slice * sliceheight + threadIdx.x % sliceheight; 21590d2215bSHong Zhang if (row < nrows) t += aval[gid] * x[acolidx[gid]]; 21690d2215bSHong Zhang if (iter == chunksperblock - 1 || (cid + 2) * BLOCKY * 32 > sliidx[start_slice + 1]) { /* last iteration or next iteration covers more than one slice */ 21790d2215bSHong Zhang int tid = threadIdx.x + threadIdx.y * 32, tidx = tid % BLOCKY, tidy = tid / BLOCKY; 21890d2215bSHong Zhang /* reduction and write to output vector */ 21990d2215bSHong Zhang #pragma unroll 22090d2215bSHong Zhang for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); } 22190d2215bSHong Zhang /* transpose layout to reduce each row using warp shfl */ 22290d2215bSHong Zhang if (threadIdx.x < sliceheight) shared[threadIdx.x * BLOCKY + threadIdx.y] = t; /* shared[threadIdx.x][threadIdx.y] = t */ 22390d2215bSHong Zhang __syncthreads(); 22490d2215bSHong Zhang if (tidy < sliceheight) t = shared[tidy * BLOCKY + tidx]; /* shared[tidy][tidx] */ 22590d2215bSHong Zhang #pragma unroll 22690d2215bSHong Zhang for (int offset = BLOCKY / 2; offset > 0; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset, BLOCKY); } 22790d2215bSHong Zhang if (tidx == 0 && tidy < sliceheight) { shared[tidy] = t; /* shared[0][tidy] = t */ } 22890d2215bSHong Zhang __syncthreads(); 22990d2215bSHong Zhang if (row < nrows && threadIdx.y == 0 && threadIdx.x < sliceheight) atomicAdd(&y[row], shared[threadIdx.x]); /* shared[0][threadIdx.x] */ 23090d2215bSHong Zhang t = 0.0; 23190d2215bSHong Zhang } 23290d2215bSHong Zhang } 23390d2215bSHong Zhang } 23490d2215bSHong Zhang } 23590d2215bSHong Zhang } 23690d2215bSHong Zhang 23790d2215bSHong Zhang /* load-balancing version. Chunksize is equal to the number of threads per block */ 23890d2215bSHong Zhang template <int BLOCKY> 23990d2215bSHong 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) 24090d2215bSHong Zhang { 24190d2215bSHong Zhang __shared__ MatScalar shared[BLOCKY * 32]; 24290d2215bSHong Zhang PetscInt gid, row, start_slice, cid; 24390d2215bSHong Zhang PetscScalar t = 0.0; 24490d2215bSHong Zhang /* copy y to z */ 24590d2215bSHong Zhang for (int iter = 0; iter < 1 + (nrows - 1) / (gridDim.x * 32 * BLOCKY); iter++) { 24690d2215bSHong Zhang gid = gridDim.x * 32 * BLOCKY * iter + blockIdx.x * BLOCKY * 32 + threadIdx.y * 32 + threadIdx.x; 24790d2215bSHong Zhang if (gid < nrows) z[gid] = y[gid]; 24890d2215bSHong Zhang } 24990d2215bSHong Zhang for (int iter = 0; iter < chunksperblock; iter++) { 25090d2215bSHong Zhang cid = blockIdx.x * chunksperblock + iter; /* chunk id */ 25190d2215bSHong Zhang if (cid < totalchunks) { 25290d2215bSHong Zhang start_slice = chunk_slice_map[cid]; /* starting slice at each iteration */ 25390d2215bSHong Zhang gid = cid * BLOCKY * 32 + threadIdx.y * 32 + threadIdx.x; 25490d2215bSHong Zhang if ((cid + 1) * BLOCKY * 32 > sliidx[start_slice + 1]) { /* this iteration covers more than one slice */ 25590d2215bSHong Zhang __shared__ PetscInt flag[BLOCKY * 2]; 25690d2215bSHong Zhang bool write; 25790d2215bSHong Zhang PetscInt slice_id = start_slice, totalslices = PetscCeilInt(nrows, sliceheight), totalentries = sliidx[totalslices]; 25890d2215bSHong Zhang /* find out the slice that this element belongs to */ 25990d2215bSHong Zhang while (gid < totalentries && gid >= sliidx[slice_id + 1]) slice_id++; 26090d2215bSHong Zhang if (threadIdx.x % 16 == 0) flag[threadIdx.y * 2 + threadIdx.x / 16] = slice_id; 26190d2215bSHong Zhang row = slice_id * sliceheight + threadIdx.x % sliceheight; 26290d2215bSHong Zhang if (row < nrows && gid < totalentries) t = aval[gid] * x[acolidx[gid]]; 26390d2215bSHong Zhang __syncthreads(); 26490d2215bSHong Zhang write = segment_scan<BLOCKY>(flag, shared, &t); 26590d2215bSHong Zhang if (row < nrows && gid < totalentries && write) atomicAdd(&z[row], t); 26690d2215bSHong Zhang t = 0.0; 26790d2215bSHong Zhang } else { /* this iteration covers only one slice */ 26890d2215bSHong Zhang row = start_slice * sliceheight + threadIdx.x % sliceheight; 26990d2215bSHong Zhang if (row < nrows) t += aval[gid] * x[acolidx[gid]]; 27090d2215bSHong Zhang if (iter == chunksperblock - 1 || (cid + 2) * BLOCKY * 32 > sliidx[start_slice + 1]) { /* last iteration or next iteration covers more than one slice */ 27190d2215bSHong Zhang int tid = threadIdx.x + threadIdx.y * 32, tidx = tid % BLOCKY, tidy = tid / BLOCKY; 27290d2215bSHong Zhang /* reduction and write to output vector */ 27390d2215bSHong Zhang #pragma unroll 27490d2215bSHong Zhang for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); } 27590d2215bSHong Zhang /* transpose layout to reduce each row using warp shfl */ 27690d2215bSHong Zhang if (threadIdx.x < sliceheight) shared[threadIdx.x * BLOCKY + threadIdx.y] = t; /* shared[threadIdx.x][threadIdx.y] = t */ 27790d2215bSHong Zhang __syncthreads(); 27890d2215bSHong Zhang if (tidy < sliceheight) t = shared[tidy * BLOCKY + tidx]; /* shared[tidy][tidx] */ 27990d2215bSHong Zhang #pragma unroll 28090d2215bSHong Zhang for (int offset = BLOCKY / 2; offset > 0; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset, BLOCKY); } 28190d2215bSHong Zhang if (tidx == 0 && tidy < sliceheight) { shared[tidy] = t; /* shared[0][tidy] = t */ } 28290d2215bSHong Zhang __syncthreads(); 28390d2215bSHong Zhang if (row < nrows && threadIdx.y == 0 && threadIdx.x < sliceheight) atomicAdd(&z[row], shared[threadIdx.x]); /* shared[0][threadIdx.x] */ 28490d2215bSHong Zhang t = 0.0; 28590d2215bSHong Zhang } 28690d2215bSHong Zhang } 28790d2215bSHong Zhang } 28890d2215bSHong Zhang } 28990d2215bSHong Zhang } 29090d2215bSHong Zhang 29107e43b41SHong Zhang /* use 1 warp per slice, suitable for small slice width */ 2924e58db63SHong 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) 2932d1451d4SHong Zhang { 29407e43b41SHong Zhang PetscInt i, row, slice_id; 29507e43b41SHong Zhang slice_id = blockIdx.x * blockDim.y + threadIdx.y; 2964e58db63SHong Zhang row = slice_id * sliceheight + threadIdx.x % sliceheight; 29707e43b41SHong Zhang double t = 0.0; 29807e43b41SHong Zhang if (row < nrows) { 29907e43b41SHong Zhang for (i = sliidx[slice_id] + threadIdx.x; i < sliidx[slice_id + 1]; i += 32) t += aval[i] * x[acolidx[i]]; 30007e43b41SHong Zhang } 3014e58db63SHong Zhang #pragma unroll 3024e58db63SHong Zhang for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); } 3034e58db63SHong Zhang if (row < nrows && threadIdx.x < sliceheight) { y[row] = t; } 30407e43b41SHong Zhang } 30507e43b41SHong Zhang 306cca9ff8bSHong Zhang /* use 1 warp per slice, suitable for small slice width */ 307cca9ff8bSHong 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) 308cca9ff8bSHong Zhang { 309cca9ff8bSHong Zhang PetscInt i, row, slice_id; 310cca9ff8bSHong Zhang slice_id = blockIdx.x * blockDim.y + threadIdx.y; 311cca9ff8bSHong Zhang row = slice_id * sliceheight + threadIdx.x % sliceheight; 312cca9ff8bSHong Zhang double t = 0.0; 313cca9ff8bSHong Zhang if (row < nrows) { 314cca9ff8bSHong Zhang for (i = sliidx[slice_id] + threadIdx.x; i < sliidx[slice_id + 1]; i += 32) t += aval[i] * x[acolidx[i]]; 315cca9ff8bSHong Zhang } 316cca9ff8bSHong Zhang #pragma unroll 317cca9ff8bSHong Zhang for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); } 318cca9ff8bSHong Zhang if (row < nrows && threadIdx.x < sliceheight) { z[row] = y[row] + t; } 319cca9ff8bSHong Zhang } 320cca9ff8bSHong Zhang 321a9dd396cSHong Zhang /*********** Kernel 2-6 are tied to slice height 16. They are kept only for performance comparison **********/ 322a9dd396cSHong Zhang 323a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel6(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 32407e43b41SHong Zhang { 32507e43b41SHong Zhang __shared__ MatScalar shared[512]; 3262d1451d4SHong Zhang PetscInt i, row, slice_id, row_in_slice; 32707e43b41SHong Zhang /* multiple threads per row. */ 3282d1451d4SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 3292d1451d4SHong Zhang if (row < nrows) { 3302d1451d4SHong Zhang slice_id = row / SLICE_HEIGHT; 3312d1451d4SHong Zhang row_in_slice = row % SLICE_HEIGHT; 3322d1451d4SHong Zhang 3332d1451d4SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 3342d1451d4SHong 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]]; 33507e43b41SHong Zhang __syncthreads(); 33607e43b41SHong Zhang if (threadIdx.y < 16) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 16) * blockDim.x + threadIdx.x]; } 33707e43b41SHong Zhang __syncthreads(); 33807e43b41SHong Zhang if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; } 3392d1451d4SHong Zhang __syncthreads(); 3402d1451d4SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 3412d1451d4SHong Zhang __syncthreads(); 3422d1451d4SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 3432d1451d4SHong Zhang __syncthreads(); 3442d1451d4SHong Zhang if (threadIdx.y < 1) { 34507e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 3462d1451d4SHong Zhang y[row] = shared[threadIdx.x]; 3472d1451d4SHong Zhang } 3482d1451d4SHong Zhang } 3492d1451d4SHong Zhang } 3502d1451d4SHong Zhang 351a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel5(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 3522d1451d4SHong Zhang { 35307e43b41SHong Zhang __shared__ MatScalar shared[512]; 3542d1451d4SHong Zhang PetscInt i, row, slice_id, row_in_slice; 35507e43b41SHong Zhang /* multiple threads per row. */ 3562d1451d4SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 3572d1451d4SHong Zhang if (row < nrows) { 3582d1451d4SHong Zhang slice_id = row / SLICE_HEIGHT; 3592d1451d4SHong Zhang row_in_slice = row % SLICE_HEIGHT; 3602d1451d4SHong Zhang 3612d1451d4SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 3622d1451d4SHong 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]]; 36307e43b41SHong Zhang __syncthreads(); 36407e43b41SHong Zhang if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; } 3652d1451d4SHong Zhang __syncthreads(); 3662d1451d4SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 3672d1451d4SHong Zhang __syncthreads(); 3682d1451d4SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 3692d1451d4SHong Zhang __syncthreads(); 3702d1451d4SHong Zhang if (threadIdx.y < 1) { 37107e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 37207e43b41SHong Zhang y[row] = shared[threadIdx.x]; 37307e43b41SHong Zhang } 37407e43b41SHong Zhang } 37507e43b41SHong Zhang } 37607e43b41SHong Zhang 377a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel4(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 37807e43b41SHong Zhang { 37907e43b41SHong Zhang __shared__ MatScalar shared[512]; 38007e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 38107e43b41SHong Zhang /* multiple threads per row. */ 38207e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 38307e43b41SHong Zhang if (row < nrows) { 38407e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 38507e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 38607e43b41SHong Zhang 38707e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 38807e43b41SHong 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]]; 38907e43b41SHong Zhang __syncthreads(); 39007e43b41SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 39107e43b41SHong Zhang __syncthreads(); 39207e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 39307e43b41SHong Zhang __syncthreads(); 39407e43b41SHong Zhang if (threadIdx.y < 1) { 39507e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 39607e43b41SHong Zhang y[row] = shared[threadIdx.x]; 39707e43b41SHong Zhang } 39807e43b41SHong Zhang } 39907e43b41SHong Zhang } 40007e43b41SHong Zhang 401a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel3(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 40207e43b41SHong Zhang { 40307e43b41SHong Zhang __shared__ MatScalar shared[512]; 40407e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 40507e43b41SHong Zhang /* multiple threads per row. */ 40607e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 40707e43b41SHong Zhang if (row < nrows) { 40807e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 40907e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 41007e43b41SHong Zhang 41107e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 41207e43b41SHong 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]]; 41307e43b41SHong Zhang __syncthreads(); 41407e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 41507e43b41SHong Zhang __syncthreads(); 41607e43b41SHong Zhang if (threadIdx.y < 1) { 41707e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 41807e43b41SHong Zhang y[row] = shared[threadIdx.x]; 41907e43b41SHong Zhang } 42007e43b41SHong Zhang } 42107e43b41SHong Zhang } 42207e43b41SHong Zhang 423a9dd396cSHong Zhang __global__ void matmult_seqsell_tiled_kernel2(PetscInt nrows, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 42407e43b41SHong Zhang { 42507e43b41SHong Zhang __shared__ MatScalar shared[512]; 42607e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 42707e43b41SHong Zhang /* multiple threads per row. */ 42807e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 42907e43b41SHong Zhang if (row < nrows) { 43007e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 43107e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 43207e43b41SHong Zhang 43307e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 43407e43b41SHong 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]]; 43507e43b41SHong Zhang __syncthreads(); 43607e43b41SHong Zhang if (threadIdx.y < 1) { 43707e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 43807e43b41SHong Zhang y[row] = shared[threadIdx.x]; 43907e43b41SHong Zhang } 44007e43b41SHong Zhang } 44107e43b41SHong Zhang } 44207e43b41SHong Zhang 443a9dd396cSHong 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) 44407e43b41SHong Zhang { 44507e43b41SHong Zhang __shared__ MatScalar shared[512]; 44607e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 44707e43b41SHong Zhang /* multiple threads per row. */ 44807e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 44907e43b41SHong Zhang if (row < nrows) { 45007e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 45107e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 45207e43b41SHong Zhang 45307e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 45407e43b41SHong 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]]; 45507e43b41SHong Zhang __syncthreads(); 45607e43b41SHong Zhang if (threadIdx.y < 16) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 16) * blockDim.x + threadIdx.x]; } 45707e43b41SHong Zhang __syncthreads(); 45807e43b41SHong Zhang if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; } 45907e43b41SHong Zhang __syncthreads(); 46007e43b41SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 46107e43b41SHong Zhang __syncthreads(); 46207e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 46307e43b41SHong Zhang __syncthreads(); 46407e43b41SHong Zhang if (threadIdx.y < 1) { 46507e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 4662d1451d4SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 4672d1451d4SHong Zhang } 4682d1451d4SHong Zhang } 4692d1451d4SHong Zhang } 47007e43b41SHong Zhang 471a9dd396cSHong 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) 47207e43b41SHong Zhang { 47307e43b41SHong Zhang __shared__ MatScalar shared[512]; 47407e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 47507e43b41SHong Zhang /* multiple threads per row. */ 47607e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 47707e43b41SHong Zhang if (row < nrows) { 47807e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 47907e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 48007e43b41SHong Zhang 48107e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 48207e43b41SHong 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]]; 48307e43b41SHong Zhang __syncthreads(); 48407e43b41SHong Zhang if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; } 48507e43b41SHong Zhang __syncthreads(); 48607e43b41SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 48707e43b41SHong Zhang __syncthreads(); 48807e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 48907e43b41SHong Zhang __syncthreads(); 49007e43b41SHong Zhang if (threadIdx.y < 1) { 49107e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 49207e43b41SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 49307e43b41SHong Zhang } 49407e43b41SHong Zhang } 49507e43b41SHong Zhang } 49607e43b41SHong Zhang 497a9dd396cSHong 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) 49807e43b41SHong Zhang { 49907e43b41SHong Zhang __shared__ MatScalar shared[512]; 50007e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 50107e43b41SHong Zhang /* multiple threads per row. */ 50207e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 50307e43b41SHong Zhang if (row < nrows) { 50407e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 50507e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 50607e43b41SHong Zhang 50707e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 50807e43b41SHong 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]]; 50907e43b41SHong Zhang __syncthreads(); 51007e43b41SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 51107e43b41SHong Zhang __syncthreads(); 51207e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 51307e43b41SHong Zhang __syncthreads(); 51407e43b41SHong Zhang if (threadIdx.y < 1) { 51507e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 51607e43b41SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 51707e43b41SHong Zhang } 51807e43b41SHong Zhang } 51907e43b41SHong Zhang } 52007e43b41SHong Zhang 521a9dd396cSHong 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) 52207e43b41SHong Zhang { 52307e43b41SHong Zhang __shared__ MatScalar shared[512]; 52407e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 52507e43b41SHong Zhang /* multiple threads per row. */ 52607e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 52707e43b41SHong Zhang if (row < nrows) { 52807e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 52907e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 53007e43b41SHong Zhang 53107e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 53207e43b41SHong 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]]; 53307e43b41SHong Zhang __syncthreads(); 53407e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 53507e43b41SHong Zhang __syncthreads(); 53607e43b41SHong Zhang if (threadIdx.y < 1) { 53707e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 53807e43b41SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 53907e43b41SHong Zhang } 54007e43b41SHong Zhang } 54107e43b41SHong Zhang } 54207e43b41SHong Zhang 543a9dd396cSHong 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) 54407e43b41SHong Zhang { 54507e43b41SHong Zhang __shared__ MatScalar shared[512]; 54607e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 54707e43b41SHong Zhang /* multiple threads per row. */ 54807e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 54907e43b41SHong Zhang if (row < nrows) { 55007e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 55107e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 55207e43b41SHong Zhang 55307e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 55407e43b41SHong 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]]; 55507e43b41SHong Zhang __syncthreads(); 55607e43b41SHong Zhang if (threadIdx.y < 1) { 55707e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 55807e43b41SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 55907e43b41SHong Zhang } 56007e43b41SHong Zhang } 5612d1451d4SHong Zhang } 5622d1451d4SHong Zhang 5632d1451d4SHong Zhang PetscErrorCode MatMult_SeqSELLCUDA(Mat A, Vec xx, Vec yy) 5642d1451d4SHong Zhang { 5652d1451d4SHong Zhang Mat_SeqSELL *a = (Mat_SeqSELL *)A->data; 5662d1451d4SHong Zhang Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr; 5672d1451d4SHong Zhang PetscScalar *y; 5682d1451d4SHong Zhang const PetscScalar *x; 569a9dd396cSHong Zhang PetscInt nrows = A->rmap->n, sliceheight = a->sliceheight; 57090d2215bSHong Zhang PetscInt chunksperblock, nchunks, *chunk_slice_map; 5712d1451d4SHong Zhang MatScalar *aval; 5722d1451d4SHong Zhang PetscInt *acolidx; 5732d1451d4SHong Zhang PetscInt *sliidx; 57407e43b41SHong Zhang PetscInt nblocks, blocksize = 512; /* blocksize must be multiple of SLICE_HEIGHT*32 */ 57507e43b41SHong Zhang dim3 block2(256, 2), block4(128, 4), block8(64, 8), block16(32, 16), block32(16, 32); 576*b921024eSHong Zhang PetscReal maxoveravg; 5772d1451d4SHong Zhang 5782d1451d4SHong Zhang PetscFunctionBegin; 5794e58db63SHong 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); 58090d2215bSHong 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); 5812d1451d4SHong Zhang PetscCall(MatSeqSELLCUDACopyToGPU(A)); 5822d1451d4SHong Zhang /* cudastruct may not be available until MatSeqSELLCUDACopyToGPU() is called */ 5832d1451d4SHong Zhang aval = cudastruct->val; 5842d1451d4SHong Zhang acolidx = cudastruct->colidx; 5852d1451d4SHong Zhang sliidx = cudastruct->sliidx; 5862d1451d4SHong Zhang 5872d1451d4SHong Zhang PetscCall(VecCUDAGetArrayRead(xx, &x)); 5882d1451d4SHong Zhang PetscCall(VecCUDAGetArrayWrite(yy, &y)); 5892d1451d4SHong Zhang PetscCall(PetscLogGpuTimeBegin()); 59007e43b41SHong Zhang 59107e43b41SHong Zhang switch (cudastruct->kernelchoice) { 59207e43b41SHong Zhang case 9: 5934e58db63SHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 5944e58db63SHong Zhang if (cudastruct->blocky == 2) { 5954e58db63SHong Zhang matmult_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 5964e58db63SHong Zhang } else if (cudastruct->blocky == 4) { 5974e58db63SHong Zhang matmult_seqsell_tiled_kernel9<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 5984e58db63SHong Zhang } else if (cudastruct->blocky == 8) { 5994e58db63SHong Zhang matmult_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 6004e58db63SHong Zhang } else if (cudastruct->blocky == 16) { 6014e58db63SHong Zhang matmult_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 6024e58db63SHong Zhang } else if (cudastruct->blocky == 32) { 6034e58db63SHong Zhang matmult_seqsell_tiled_kernel9<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 604cca9ff8bSHong Zhang } else { 605cca9ff8bSHong Zhang matmult_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 60607e43b41SHong Zhang } 60707e43b41SHong Zhang break; 60807e43b41SHong Zhang case 7: 6094e58db63SHong Zhang nblocks = 1 + (nrows - 1) / (2 * sliceheight); 6104e58db63SHong Zhang if (cudastruct->blocky == 2) { 6114e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 6124e58db63SHong Zhang } else if (cudastruct->blocky == 4) { 6134e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 6144e58db63SHong Zhang } else if (cudastruct->blocky == 8) { 6154e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 6164e58db63SHong Zhang } else if (cudastruct->blocky == 16) { 6174e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 6184e58db63SHong Zhang } else if (cudastruct->blocky == 32) { 6194e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 6204e58db63SHong Zhang } else { 6214e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 6224e58db63SHong Zhang } 62307e43b41SHong Zhang break; 62407e43b41SHong Zhang case 6: 62507e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 32); /* 1 slice per block if blocksize=512 */ 626a9dd396cSHong Zhang matmult_seqsell_tiled_kernel6<<<nblocks, block32>>>(nrows, acolidx, aval, sliidx, x, y); 62707e43b41SHong Zhang break; 62807e43b41SHong Zhang case 5: 62907e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 16); /* 2 slices per block if blocksize=512*/ 630a9dd396cSHong Zhang matmult_seqsell_tiled_kernel5<<<nblocks, block16>>>(nrows, acolidx, aval, sliidx, x, y); 63107e43b41SHong Zhang break; 63207e43b41SHong Zhang case 4: 63307e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 8); /* 4 slices per block if blocksize=512 */ 634a9dd396cSHong Zhang matmult_seqsell_tiled_kernel4<<<nblocks, block8>>>(nrows, acolidx, aval, sliidx, x, y); 63507e43b41SHong Zhang break; 63607e43b41SHong Zhang case 3: 63707e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 4); /* 8 slices per block if blocksize=512 */ 638a9dd396cSHong Zhang matmult_seqsell_tiled_kernel3<<<nblocks, block4>>>(nrows, acolidx, aval, sliidx, x, y); 63907e43b41SHong Zhang break; 64007e43b41SHong Zhang case 2: /* 16 slices per block if blocksize=512 */ 64107e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 2); 642a9dd396cSHong Zhang matmult_seqsell_tiled_kernel2<<<nblocks, block2>>>(nrows, acolidx, aval, sliidx, x, y); 64307e43b41SHong Zhang break; 64407e43b41SHong Zhang case 1: /* 32 slices per block if blocksize=512 */ 64507e43b41SHong Zhang nblocks = 1 + (nrows - 1) / blocksize; 6464e58db63SHong Zhang matmult_seqsell_basic_kernel<<<nblocks, blocksize>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 64707e43b41SHong Zhang break; 64807e43b41SHong Zhang case 0: 64990d2215bSHong Zhang maxoveravg = a->maxslicewidth / a->avgslicewidth; 65090d2215bSHong Zhang if (maxoveravg > 12.0 && maxoveravg / nrows > 0.001) { /* important threshold */ 65190d2215bSHong Zhang /* each block handles approximately one slice */ 65290d2215bSHong Zhang PetscInt blocky = a->chunksize / 32; 65390d2215bSHong Zhang nchunks = cudastruct->totalchunks; 65490d2215bSHong Zhang chunksperblock = cudastruct->chunksperblock ? cudastruct->chunksperblock : 1 + (cudastruct->totalentries / cudastruct->totalslices - 1) / a->chunksize; 65590d2215bSHong Zhang nblocks = 1 + (nchunks - 1) / chunksperblock; 65690d2215bSHong Zhang chunk_slice_map = cudastruct->chunk_slice_map; 65790d2215bSHong Zhang if (blocky == 2) { 65890d2215bSHong Zhang matmult_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y); 65990d2215bSHong Zhang } else if (blocky == 4) { 66090d2215bSHong Zhang matmult_seqsell_tiled_kernel8<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y); 66190d2215bSHong Zhang } else if (blocky == 8) { 66290d2215bSHong Zhang matmult_seqsell_tiled_kernel8<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y); 66390d2215bSHong Zhang } else if (blocky == 16) { 66490d2215bSHong Zhang matmult_seqsell_tiled_kernel8<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y); 66590d2215bSHong Zhang } else if (blocky == 32) { 66690d2215bSHong Zhang matmult_seqsell_tiled_kernel8<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y); 66790d2215bSHong Zhang } else { 66890d2215bSHong Zhang matmult_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y); 66990d2215bSHong Zhang } 6702d1451d4SHong Zhang } else { 671cca9ff8bSHong Zhang PetscInt avgslicesize = sliceheight * a->avgslicewidth; 67290d2215bSHong Zhang if (avgslicesize <= 432) { 67390d2215bSHong Zhang if (sliceheight * a->maxslicewidth < 2048 && nrows > 100000) { 674cca9ff8bSHong Zhang nblocks = 1 + (nrows - 1) / (2 * sliceheight); /* two slices per block */ 675cca9ff8bSHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 67690d2215bSHong Zhang } else { 677cca9ff8bSHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 678cca9ff8bSHong Zhang matmult_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 67990d2215bSHong Zhang } 680cca9ff8bSHong Zhang } else if (avgslicesize <= 2400) { 681cca9ff8bSHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 682cca9ff8bSHong Zhang matmult_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 683cca9ff8bSHong Zhang } else { 6844e58db63SHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 6854e58db63SHong Zhang matmult_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 6862d1451d4SHong Zhang } 6872d1451d4SHong Zhang } 688cca9ff8bSHong Zhang break; 689cca9ff8bSHong Zhang } 6902d1451d4SHong Zhang PetscCall(PetscLogGpuTimeEnd()); 6912d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayRead(xx, &x)); 6922d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayWrite(yy, &y)); 6932d1451d4SHong Zhang PetscCall(PetscLogGpuFlops(2.0 * a->nz - a->nonzerorowcnt)); 6942d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 6952d1451d4SHong Zhang } 6962d1451d4SHong Zhang 6972d1451d4SHong Zhang PetscErrorCode MatMultAdd_SeqSELLCUDA(Mat A, Vec xx, Vec yy, Vec zz) 6982d1451d4SHong Zhang { 6992d1451d4SHong Zhang Mat_SeqSELL *a = (Mat_SeqSELL *)A->data; 7002d1451d4SHong Zhang Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr; 7012d1451d4SHong Zhang PetscScalar *z; 7022d1451d4SHong Zhang const PetscScalar *y, *x; 703a9dd396cSHong Zhang PetscInt nrows = A->rmap->n, sliceheight = a->sliceheight; 70490d2215bSHong Zhang PetscInt chunksperblock, nchunks, *chunk_slice_map; 7052d1451d4SHong Zhang MatScalar *aval = cudastruct->val; 7062d1451d4SHong Zhang PetscInt *acolidx = cudastruct->colidx; 7072d1451d4SHong Zhang PetscInt *sliidx = cudastruct->sliidx; 70890d2215bSHong Zhang PetscReal maxoveravg; 7092d1451d4SHong Zhang 7102d1451d4SHong Zhang PetscFunctionBegin; 71190d2215bSHong 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); 71290d2215bSHong 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); 7132d1451d4SHong Zhang PetscCall(MatSeqSELLCUDACopyToGPU(A)); 7142d1451d4SHong Zhang if (a->nz) { 71590d2215bSHong Zhang PetscInt blocky = cudastruct->blocky, nblocks, blocksize = 512; 71607e43b41SHong Zhang dim3 block2(256, 2), block4(128, 4), block8(64, 8), block16(32, 16), block32(16, 32); 7172d1451d4SHong Zhang PetscCall(VecCUDAGetArrayRead(xx, &x)); 7182d1451d4SHong Zhang PetscCall(VecCUDAGetArrayRead(yy, &y)); 7192d1451d4SHong Zhang PetscCall(VecCUDAGetArrayWrite(zz, &z)); 7202d1451d4SHong Zhang PetscCall(PetscLogGpuTimeBegin()); 72107e43b41SHong Zhang 72207e43b41SHong Zhang switch (cudastruct->kernelchoice) { 723cca9ff8bSHong Zhang case 9: 724cca9ff8bSHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 72590d2215bSHong Zhang if (blocky == 2) { 726cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 72790d2215bSHong Zhang } else if (blocky == 4) { 728cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel9<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 72990d2215bSHong Zhang } else if (blocky == 8) { 730cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 73190d2215bSHong Zhang } else if (blocky == 16) { 732cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 73390d2215bSHong Zhang } else if (blocky == 32) { 734cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel9<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 735cca9ff8bSHong Zhang } else { 736cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 737cca9ff8bSHong Zhang } 738cca9ff8bSHong Zhang break; 73990d2215bSHong Zhang case 8: 74090d2215bSHong Zhang /* each block handles approximately one slice */ 74190d2215bSHong Zhang nchunks = cudastruct->totalchunks; 74290d2215bSHong Zhang blocky = a->chunksize / 32; 74390d2215bSHong Zhang chunksperblock = cudastruct->chunksperblock ? cudastruct->chunksperblock : 1 + (cudastruct->totalentries / cudastruct->totalslices - 1) / a->chunksize; 74490d2215bSHong Zhang nblocks = 1 + (nchunks - 1) / chunksperblock; 74590d2215bSHong Zhang chunk_slice_map = cudastruct->chunk_slice_map; 74690d2215bSHong Zhang if (blocky == 2) { 74790d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 74890d2215bSHong Zhang } else if (blocky == 4) { 74990d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 75090d2215bSHong Zhang } else if (blocky == 8) { 75190d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 75290d2215bSHong Zhang } else if (blocky == 16) { 75390d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 75490d2215bSHong Zhang } else if (blocky == 32) { 75590d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 75690d2215bSHong Zhang } else { 75790d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 75890d2215bSHong Zhang } 75990d2215bSHong Zhang break; 760cca9ff8bSHong Zhang case 7: 761cca9ff8bSHong Zhang nblocks = 1 + (nrows - 1) / (2 * sliceheight); 76290d2215bSHong Zhang if (blocky == 2) { 763cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 76490d2215bSHong Zhang } else if (blocky == 4) { 765cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 76690d2215bSHong Zhang } else if (blocky == 8) { 767cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 76890d2215bSHong Zhang } else if (blocky == 16) { 769cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 77090d2215bSHong Zhang } else if (blocky == 32) { 771cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 772cca9ff8bSHong Zhang } else { 773cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 774cca9ff8bSHong Zhang } 775cca9ff8bSHong Zhang break; 77607e43b41SHong Zhang case 6: 77707e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 32); 778a9dd396cSHong Zhang matmultadd_seqsell_tiled_kernel6<<<nblocks, block32>>>(nrows, acolidx, aval, sliidx, x, y, z); 77907e43b41SHong Zhang break; 78007e43b41SHong Zhang case 5: 78107e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 16); 782a9dd396cSHong Zhang matmultadd_seqsell_tiled_kernel5<<<nblocks, block16>>>(nrows, acolidx, aval, sliidx, x, y, z); 78307e43b41SHong Zhang break; 78407e43b41SHong Zhang case 4: 78507e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 8); 786a9dd396cSHong Zhang matmultadd_seqsell_tiled_kernel4<<<nblocks, block8>>>(nrows, acolidx, aval, sliidx, x, y, z); 78707e43b41SHong Zhang break; 78807e43b41SHong Zhang case 3: 78907e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 4); 790a9dd396cSHong Zhang matmultadd_seqsell_tiled_kernel3<<<nblocks, block4>>>(nrows, acolidx, aval, sliidx, x, y, z); 79107e43b41SHong Zhang break; 79207e43b41SHong Zhang case 2: 79307e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 2); 794a9dd396cSHong Zhang matmultadd_seqsell_tiled_kernel2<<<nblocks, block2>>>(nrows, acolidx, aval, sliidx, x, y, z); 79507e43b41SHong Zhang break; 79607e43b41SHong Zhang case 1: 79707e43b41SHong Zhang nblocks = 1 + (nrows - 1) / blocksize; 7984e58db63SHong Zhang matmultadd_seqsell_basic_kernel<<<nblocks, blocksize>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 79907e43b41SHong Zhang break; 800cca9ff8bSHong Zhang case 0: 80190d2215bSHong Zhang maxoveravg = a->maxslicewidth / a->avgslicewidth; 80290d2215bSHong Zhang if (maxoveravg > 12.0 && maxoveravg / nrows > 0.001) { /* important threshold */ 80390d2215bSHong Zhang /* each block handles approximately one slice */ 80490d2215bSHong Zhang nchunks = cudastruct->totalchunks; 80590d2215bSHong Zhang blocky = a->chunksize / 32; 80690d2215bSHong Zhang chunksperblock = cudastruct->chunksperblock ? cudastruct->chunksperblock : 1 + (cudastruct->totalentries / cudastruct->totalslices - 1) / a->chunksize; 80790d2215bSHong Zhang nblocks = 1 + (nchunks - 1) / chunksperblock; 80890d2215bSHong Zhang chunk_slice_map = cudastruct->chunk_slice_map; 80990d2215bSHong Zhang if (blocky == 2) { 81090d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 81190d2215bSHong Zhang } else if (blocky == 4) { 81290d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 81390d2215bSHong Zhang } else if (blocky == 8) { 81490d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 81590d2215bSHong Zhang } else if (blocky == 16) { 81690d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 81790d2215bSHong Zhang } else if (blocky == 32) { 81890d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 81990d2215bSHong Zhang } else { 82090d2215bSHong Zhang matmultadd_seqsell_tiled_kernel8<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, chunksperblock, nchunks, chunk_slice_map, acolidx, aval, sliidx, x, y, z); 82190d2215bSHong Zhang } 822cca9ff8bSHong Zhang } else { 823cca9ff8bSHong Zhang PetscInt avgslicesize = sliceheight * a->avgslicewidth; 82490d2215bSHong Zhang if (avgslicesize <= 432) { 82590d2215bSHong Zhang if (sliceheight * a->maxslicewidth < 2048 && nrows > 100000) { 826cca9ff8bSHong Zhang nblocks = 1 + (nrows - 1) / (2 * sliceheight); /* two slices per block */ 827cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 82890d2215bSHong Zhang } else { 829cca9ff8bSHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 830cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 83190d2215bSHong Zhang } 832cca9ff8bSHong Zhang } else if (avgslicesize <= 2400) { 833cca9ff8bSHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 834cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 835cca9ff8bSHong Zhang } else { 836cca9ff8bSHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 837cca9ff8bSHong Zhang matmultadd_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 838cca9ff8bSHong Zhang } 839cca9ff8bSHong Zhang } 84007e43b41SHong Zhang break; 8412d1451d4SHong Zhang } 8422d1451d4SHong Zhang PetscCall(PetscLogGpuTimeEnd()); 8432d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayRead(xx, &x)); 8442d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayRead(yy, &y)); 8452d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayWrite(zz, &z)); 8462d1451d4SHong Zhang PetscCall(PetscLogGpuFlops(2.0 * a->nz)); 8472d1451d4SHong Zhang } else { 8482d1451d4SHong Zhang PetscCall(VecCopy(yy, zz)); 8492d1451d4SHong Zhang } 8502d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 8512d1451d4SHong Zhang } 8522d1451d4SHong Zhang 8532d1451d4SHong Zhang static PetscErrorCode MatSetFromOptions_SeqSELLCUDA(Mat A, PetscOptionItems *PetscOptionsObject) 8542d1451d4SHong Zhang { 85507e43b41SHong Zhang Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr; 8564e58db63SHong Zhang PetscInt kernel, blocky; 85707e43b41SHong Zhang PetscBool flg; 85807e43b41SHong Zhang 8592d1451d4SHong Zhang PetscFunctionBegin; 8602d1451d4SHong Zhang PetscOptionsHeadBegin(PetscOptionsObject, "SeqSELLCUDA options"); 86190d2215bSHong Zhang PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_sell_spmv_cuda_blocky", &blocky, &flg)); 86290d2215bSHong Zhang if (flg) { 86390d2215bSHong 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); 86490d2215bSHong Zhang cudastruct->blocky = blocky; 86590d2215bSHong Zhang } 86607e43b41SHong Zhang PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_sell_spmv_cuda_kernel", &kernel, &flg)); 86707e43b41SHong Zhang if (flg) { 86807e43b41SHong 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); 86907e43b41SHong Zhang cudastruct->kernelchoice = kernel; 87090d2215bSHong Zhang if (kernel == 8) { PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_sell_spmv_cuda_chunksperblock", &cudastruct->chunksperblock, &flg)); } 8714e58db63SHong Zhang } 8722d1451d4SHong Zhang PetscOptionsHeadEnd(); 8732d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 8742d1451d4SHong Zhang } 8752d1451d4SHong Zhang 87607e43b41SHong Zhang PETSC_INTERN PetscErrorCode MatAssemblyEnd_SpMV_Preprocessing_Private(Mat A) 87707e43b41SHong Zhang { 87807e43b41SHong Zhang Mat_SeqSELL *a = (Mat_SeqSELL *)A->data; 87907e43b41SHong Zhang 88090d2215bSHong Zhang PetscFunctionBegin; 88107e43b41SHong Zhang PetscCall(MatSeqSELLGetAvgSliceWidth(A, &a->avgslicewidth)); 88207e43b41SHong Zhang PetscCall(MatSeqSELLGetMaxSliceWidth(A, &a->maxslicewidth)); 88307e43b41SHong Zhang PetscCall(MatSeqSELLGetFillRatio(A, &a->fillratio)); 88407e43b41SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 88507e43b41SHong Zhang } 88607e43b41SHong Zhang 8872d1451d4SHong Zhang static PetscErrorCode MatAssemblyEnd_SeqSELLCUDA(Mat A, MatAssemblyType mode) 8882d1451d4SHong Zhang { 8892d1451d4SHong Zhang PetscFunctionBegin; 8902d1451d4SHong Zhang PetscCall(MatAssemblyEnd_SeqSELL(A, mode)); 89107e43b41SHong Zhang PetscCall(MatAssemblyEnd_SpMV_Preprocessing_Private(A)); 8922d1451d4SHong Zhang if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS); 8932d1451d4SHong Zhang if (A->factortype == MAT_FACTOR_NONE) { PetscCall(MatSeqSELLCUDACopyToGPU(A)); } 8942d1451d4SHong Zhang A->ops->mult = MatMult_SeqSELLCUDA; 8952d1451d4SHong Zhang A->ops->multadd = MatMultAdd_SeqSELLCUDA; 8962d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 8972d1451d4SHong Zhang } 8982d1451d4SHong Zhang 8992d1451d4SHong Zhang static PetscErrorCode MatDestroy_SeqSELLCUDA(Mat A) 9002d1451d4SHong Zhang { 9012d1451d4SHong Zhang PetscFunctionBegin; 9022d1451d4SHong Zhang if (A->factortype == MAT_FACTOR_NONE) { 9032d1451d4SHong Zhang if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) { PetscCall(MatSeqSELLCUDA_Destroy((Mat_SeqSELLCUDA **)&A->spptr)); } 9042d1451d4SHong Zhang } 9052d1451d4SHong Zhang PetscCall(MatDestroy_SeqSELL(A)); 9062d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 9072d1451d4SHong Zhang } 9082d1451d4SHong Zhang 9092d1451d4SHong Zhang PETSC_INTERN PetscErrorCode MatConvert_SeqSELL_SeqSELLCUDA(Mat); 9102d1451d4SHong Zhang static PetscErrorCode MatDuplicate_SeqSELLCUDA(Mat A, MatDuplicateOption cpvalues, Mat *B) 9112d1451d4SHong Zhang { 9122d1451d4SHong Zhang PetscFunctionBegin; 9132d1451d4SHong Zhang PetscCall(MatDuplicate_SeqSELL(A, cpvalues, B)); 9142d1451d4SHong Zhang PetscCall(MatConvert_SeqSELL_SeqSELLCUDA(*B)); 9152d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 9162d1451d4SHong Zhang } 9172d1451d4SHong Zhang 9182d1451d4SHong Zhang PETSC_EXTERN PetscErrorCode MatConvert_SeqSELL_SeqSELLCUDA(Mat B) 9192d1451d4SHong Zhang { 9202d1451d4SHong Zhang Mat_SeqSELLCUDA *cudastruct; 9212d1451d4SHong Zhang 9222d1451d4SHong Zhang PetscFunctionBegin; 9232d1451d4SHong Zhang PetscCall(PetscFree(B->defaultvectype)); 9242d1451d4SHong Zhang PetscCall(PetscStrallocpy(VECCUDA, &B->defaultvectype)); 9252d1451d4SHong Zhang 9262d1451d4SHong Zhang if (!B->spptr) { 9272d1451d4SHong Zhang if (B->factortype == MAT_FACTOR_NONE) { 9282d1451d4SHong Zhang PetscCall(PetscNew(&cudastruct)); 9292d1451d4SHong Zhang B->spptr = cudastruct; 9302d1451d4SHong Zhang } 9312d1451d4SHong Zhang } 9322d1451d4SHong Zhang 9332d1451d4SHong Zhang B->ops->assemblyend = MatAssemblyEnd_SeqSELLCUDA; 9342d1451d4SHong Zhang B->ops->destroy = MatDestroy_SeqSELLCUDA; 9352d1451d4SHong Zhang B->ops->setfromoptions = MatSetFromOptions_SeqSELLCUDA; 9362d1451d4SHong Zhang B->ops->mult = MatMult_SeqSELLCUDA; 9372d1451d4SHong Zhang B->ops->multadd = MatMultAdd_SeqSELLCUDA; 9382d1451d4SHong Zhang B->ops->duplicate = MatDuplicate_SeqSELLCUDA; 9392d1451d4SHong Zhang 94007e43b41SHong Zhang /* No need to assemble SeqSELL, but need to do the preprocessing for SpMV */ 94107e43b41SHong Zhang PetscCall(MatAssemblyEnd_SpMV_Preprocessing_Private(B)); 94207e43b41SHong Zhang 9432d1451d4SHong Zhang PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQSELLCUDA)); 9442d1451d4SHong Zhang B->offloadmask = PETSC_OFFLOAD_UNALLOCATED; 9452d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 9462d1451d4SHong Zhang } 9472d1451d4SHong Zhang 9482d1451d4SHong Zhang PETSC_EXTERN PetscErrorCode MatCreate_SeqSELLCUDA(Mat B) 9492d1451d4SHong Zhang { 9502d1451d4SHong Zhang PetscFunctionBegin; 9512d1451d4SHong Zhang PetscCall(MatCreate_SeqSELL(B)); 9522d1451d4SHong Zhang PetscCall(MatConvert_SeqSELL_SeqSELLCUDA(B)); 9532d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 9542d1451d4SHong Zhang } 955