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 { 92d1451d4SHong Zhang PetscInt *colidx; /* column index */ 102d1451d4SHong Zhang MatScalar *val; 112d1451d4SHong Zhang PetscInt *sliidx; 122d1451d4SHong Zhang PetscInt nonzerostate; 1307e43b41SHong Zhang PetscInt kernelchoice; 14*4e58db63SHong Zhang PetscInt blocky; 152d1451d4SHong Zhang } Mat_SeqSELLCUDA; 162d1451d4SHong Zhang 172d1451d4SHong Zhang static PetscErrorCode MatSeqSELLCUDA_Destroy(Mat_SeqSELLCUDA **cudastruct) 182d1451d4SHong Zhang { 192d1451d4SHong Zhang PetscFunctionBegin; 202d1451d4SHong Zhang if (*cudastruct) { 212d1451d4SHong Zhang if ((*cudastruct)->colidx) { PetscCallCUDA(cudaFree((*cudastruct)->colidx)); } 222d1451d4SHong Zhang if ((*cudastruct)->val) { PetscCallCUDA(cudaFree((*cudastruct)->val)); } 232d1451d4SHong Zhang if ((*cudastruct)->sliidx) { PetscCallCUDA(cudaFree((*cudastruct)->sliidx)); } 242d1451d4SHong Zhang PetscCall(PetscFree(*cudastruct)); 252d1451d4SHong Zhang } 262d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 272d1451d4SHong Zhang } 282d1451d4SHong Zhang 292d1451d4SHong Zhang static PetscErrorCode MatSeqSELLCUDACopyToGPU(Mat A) 302d1451d4SHong Zhang { 312d1451d4SHong Zhang Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr; 322d1451d4SHong Zhang Mat_SeqSELL *a = (Mat_SeqSELL *)A->data; 332d1451d4SHong Zhang 342d1451d4SHong Zhang PetscFunctionBegin; 352d1451d4SHong Zhang if (A->offloadmask == PETSC_OFFLOAD_UNALLOCATED || A->offloadmask == PETSC_OFFLOAD_CPU) { 362d1451d4SHong Zhang PetscCall(PetscLogEventBegin(MAT_CUDACopyToGPU, A, 0, 0, 0)); 372d1451d4SHong Zhang if (A->assembled && A->nonzerostate == cudastruct->nonzerostate) { 382d1451d4SHong Zhang /* copy values only */ 392d1451d4SHong Zhang PetscCallCUDA(cudaMemcpy(cudastruct->val, a->val, a->sliidx[a->totalslices] * sizeof(MatScalar), cudaMemcpyHostToDevice)); 402d1451d4SHong Zhang PetscCall(PetscLogCpuToGpu(a->sliidx[a->totalslices] * (sizeof(MatScalar)))); 412d1451d4SHong Zhang } else { 422d1451d4SHong Zhang if (cudastruct->colidx) { PetscCallCUDA(cudaFree(cudastruct->colidx)); } 432d1451d4SHong Zhang if (cudastruct->val) { PetscCallCUDA(cudaFree(cudastruct->val)); } 442d1451d4SHong Zhang if (cudastruct->sliidx) { PetscCallCUDA(cudaFree(cudastruct->sliidx)); } 452d1451d4SHong Zhang PetscCallCUDA(cudaMalloc((void **)&(cudastruct->colidx), a->maxallocmat * sizeof(PetscInt))); 462d1451d4SHong Zhang PetscCallCUDA(cudaMalloc((void **)&(cudastruct->val), a->maxallocmat * sizeof(MatScalar))); 472d1451d4SHong Zhang /* copy values, nz or maxallocmat? */ 482d1451d4SHong Zhang PetscCallCUDA(cudaMemcpy(cudastruct->colidx, a->colidx, a->sliidx[a->totalslices] * sizeof(PetscInt), cudaMemcpyHostToDevice)); 492d1451d4SHong Zhang PetscCallCUDA(cudaMemcpy(cudastruct->val, a->val, a->sliidx[a->totalslices] * sizeof(MatScalar), cudaMemcpyHostToDevice)); 5007e43b41SHong Zhang 5107e43b41SHong Zhang PetscCallCUDA(cudaMalloc((void **)&(cudastruct->sliidx), (a->totalslices + 1) * sizeof(PetscInt))); 5207e43b41SHong Zhang PetscCallCUDA(cudaMemcpy(cudastruct->sliidx, a->sliidx, (a->totalslices + 1) * sizeof(PetscInt), cudaMemcpyHostToDevice)); 532d1451d4SHong Zhang PetscCall(PetscLogCpuToGpu(a->sliidx[a->totalslices] * (sizeof(MatScalar) + sizeof(PetscInt)) + (a->totalslices + 1) * sizeof(PetscInt))); 542d1451d4SHong Zhang cudastruct->nonzerostate = A->nonzerostate; 552d1451d4SHong Zhang } 562d1451d4SHong Zhang PetscCallCUDA(WaitForCUDA()); 572d1451d4SHong Zhang PetscCall(PetscLogEventEnd(MAT_CUDACopyToGPU, A, 0, 0, 0)); 582d1451d4SHong Zhang A->offloadmask = PETSC_OFFLOAD_BOTH; 592d1451d4SHong Zhang } 602d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 612d1451d4SHong Zhang } 622d1451d4SHong Zhang 63*4e58db63SHong 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) 642d1451d4SHong Zhang { 652d1451d4SHong Zhang PetscInt i, row, slice_id, row_in_slice; 662d1451d4SHong Zhang MatScalar sum; 672d1451d4SHong Zhang /* one thread per row. */ 682d1451d4SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 692d1451d4SHong Zhang if (row < nrows) { 70*4e58db63SHong Zhang slice_id = row / sliceheight; 71*4e58db63SHong Zhang row_in_slice = row % sliceheight; 722d1451d4SHong Zhang sum = 0.0; 73*4e58db63SHong Zhang for (i = sliidx[slice_id] + row_in_slice; i < sliidx[slice_id + 1]; i += sliceheight) sum += aval[i] * x[acolidx[i]]; 742d1451d4SHong Zhang y[row] = sum; 752d1451d4SHong Zhang } 762d1451d4SHong Zhang } 772d1451d4SHong Zhang 78*4e58db63SHong 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) 792d1451d4SHong Zhang { 802d1451d4SHong Zhang PetscInt i, row, slice_id, row_in_slice; 812d1451d4SHong Zhang MatScalar sum; 822d1451d4SHong Zhang /* one thread per row. */ 832d1451d4SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 842d1451d4SHong Zhang if (row < nrows) { 85*4e58db63SHong Zhang slice_id = row / sliceheight; 86*4e58db63SHong Zhang row_in_slice = row % sliceheight; 872d1451d4SHong Zhang sum = 0.0; 88*4e58db63SHong Zhang for (i = sliidx[slice_id] + row_in_slice; i < sliidx[slice_id + 1]; i += sliceheight) sum += aval[i] * x[acolidx[i]]; 892d1451d4SHong Zhang z[row] = y[row] + sum; 902d1451d4SHong Zhang } 912d1451d4SHong Zhang } 9207e43b41SHong Zhang 9307e43b41SHong Zhang /* use 1 block per slice, suitable for large slice width */ 9407e43b41SHong Zhang template <int BLOCKY> 95*4e58db63SHong 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) 9607e43b41SHong Zhang { 97*4e58db63SHong Zhang __shared__ MatScalar shared[32][BLOCKY]; 98*4e58db63SHong Zhang PetscInt i, row, slice_id = blockIdx.x; 99*4e58db63SHong Zhang int tid = threadIdx.x + threadIdx.y * 32; 10007e43b41SHong Zhang /* transposed index */ 10107e43b41SHong Zhang int tidx = tid % BLOCKY; 10207e43b41SHong Zhang int tidy = tid / BLOCKY; 10307e43b41SHong Zhang PetscScalar t = 0.0; 104*4e58db63SHong Zhang 105*4e58db63SHong Zhang row = slice_id * sliceheight + threadIdx.x % sliceheight; 10607e43b41SHong Zhang if (row < nrows) { 107*4e58db63SHong 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]]; 1082d1451d4SHong Zhang } 109*4e58db63SHong Zhang #pragma unroll 110*4e58db63SHong Zhang for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); } 11107e43b41SHong Zhang /* transpose layout to reduce each row using warp shfl */ 11207e43b41SHong Zhang shared[threadIdx.x][threadIdx.y] = t; 11307e43b41SHong Zhang __syncthreads(); 11407e43b41SHong Zhang t = shared[tidy][tidx]; 11507e43b41SHong Zhang #pragma unroll 11607e43b41SHong Zhang for (int offset = BLOCKY / 2; offset > 0; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset, BLOCKY); } 117*4e58db63SHong Zhang if (tidx == 0 && tidy < sliceheight) { shared[0][tidy] = t; } 11807e43b41SHong Zhang __syncthreads(); 119*4e58db63SHong Zhang if (row < nrows && threadIdx.y == 0 && threadIdx.x < sliceheight) { y[row] = shared[0][threadIdx.x]; } 1202d1451d4SHong Zhang } 1212d1451d4SHong Zhang 12207e43b41SHong Zhang /* use 1 warp per slice, suitable for small slice width */ 123*4e58db63SHong 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) 1242d1451d4SHong Zhang { 12507e43b41SHong Zhang PetscInt i, row, slice_id; 12607e43b41SHong Zhang slice_id = blockIdx.x * blockDim.y + threadIdx.y; 127*4e58db63SHong Zhang row = slice_id * sliceheight + threadIdx.x % sliceheight; 12807e43b41SHong Zhang double t = 0.0; 12907e43b41SHong Zhang if (row < nrows) { 13007e43b41SHong Zhang for (i = sliidx[slice_id] + threadIdx.x; i < sliidx[slice_id + 1]; i += 32) t += aval[i] * x[acolidx[i]]; 13107e43b41SHong Zhang } 132*4e58db63SHong Zhang #pragma unroll 133*4e58db63SHong Zhang for (int offset = 16; offset >= sliceheight; offset /= 2) { t += __shfl_down_sync(0xffffffff, t, offset); } 134*4e58db63SHong Zhang if (row < nrows && threadIdx.x < sliceheight) { y[row] = t; } 13507e43b41SHong Zhang } 13607e43b41SHong Zhang 13707e43b41SHong Zhang __global__ void matmult_seqsell_tiled_kernel6(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 13807e43b41SHong Zhang { 13907e43b41SHong Zhang __shared__ MatScalar shared[512]; 1402d1451d4SHong Zhang PetscInt i, row, slice_id, row_in_slice; 14107e43b41SHong Zhang /* multiple threads per row. */ 1422d1451d4SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 1432d1451d4SHong Zhang if (row < nrows) { 1442d1451d4SHong Zhang slice_id = row / SLICE_HEIGHT; 1452d1451d4SHong Zhang row_in_slice = row % SLICE_HEIGHT; 1462d1451d4SHong Zhang 1472d1451d4SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 1482d1451d4SHong 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]]; 14907e43b41SHong Zhang __syncthreads(); 15007e43b41SHong Zhang if (threadIdx.y < 16) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 16) * blockDim.x + threadIdx.x]; } 15107e43b41SHong Zhang __syncthreads(); 15207e43b41SHong Zhang if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; } 1532d1451d4SHong Zhang __syncthreads(); 1542d1451d4SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 1552d1451d4SHong Zhang __syncthreads(); 1562d1451d4SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 1572d1451d4SHong Zhang __syncthreads(); 1582d1451d4SHong Zhang if (threadIdx.y < 1) { 15907e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 1602d1451d4SHong Zhang y[row] = shared[threadIdx.x]; 1612d1451d4SHong Zhang } 1622d1451d4SHong Zhang } 1632d1451d4SHong Zhang } 1642d1451d4SHong Zhang 16507e43b41SHong Zhang __global__ void matmult_seqsell_tiled_kernel5(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 1662d1451d4SHong Zhang { 16707e43b41SHong Zhang __shared__ MatScalar shared[512]; 1682d1451d4SHong Zhang PetscInt i, row, slice_id, row_in_slice; 16907e43b41SHong Zhang /* multiple threads per row. */ 1702d1451d4SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 1712d1451d4SHong Zhang if (row < nrows) { 1722d1451d4SHong Zhang slice_id = row / SLICE_HEIGHT; 1732d1451d4SHong Zhang row_in_slice = row % SLICE_HEIGHT; 1742d1451d4SHong Zhang 1752d1451d4SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 1762d1451d4SHong 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]]; 17707e43b41SHong Zhang __syncthreads(); 17807e43b41SHong Zhang if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; } 1792d1451d4SHong Zhang __syncthreads(); 1802d1451d4SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 1812d1451d4SHong Zhang __syncthreads(); 1822d1451d4SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 1832d1451d4SHong Zhang __syncthreads(); 1842d1451d4SHong Zhang if (threadIdx.y < 1) { 18507e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 18607e43b41SHong Zhang y[row] = shared[threadIdx.x]; 18707e43b41SHong Zhang } 18807e43b41SHong Zhang } 18907e43b41SHong Zhang } 19007e43b41SHong Zhang 19107e43b41SHong Zhang __global__ void matmult_seqsell_tiled_kernel4(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 19207e43b41SHong Zhang { 19307e43b41SHong Zhang __shared__ MatScalar shared[512]; 19407e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 19507e43b41SHong Zhang /* multiple threads per row. */ 19607e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 19707e43b41SHong Zhang if (row < nrows) { 19807e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 19907e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 20007e43b41SHong Zhang 20107e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 20207e43b41SHong 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]]; 20307e43b41SHong Zhang __syncthreads(); 20407e43b41SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 20507e43b41SHong Zhang __syncthreads(); 20607e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 20707e43b41SHong Zhang __syncthreads(); 20807e43b41SHong Zhang if (threadIdx.y < 1) { 20907e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 21007e43b41SHong Zhang y[row] = shared[threadIdx.x]; 21107e43b41SHong Zhang } 21207e43b41SHong Zhang } 21307e43b41SHong Zhang } 21407e43b41SHong Zhang 21507e43b41SHong Zhang __global__ void matmult_seqsell_tiled_kernel3(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 21607e43b41SHong Zhang { 21707e43b41SHong Zhang __shared__ MatScalar shared[512]; 21807e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 21907e43b41SHong Zhang /* multiple threads per row. */ 22007e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 22107e43b41SHong Zhang if (row < nrows) { 22207e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 22307e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 22407e43b41SHong Zhang 22507e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 22607e43b41SHong 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]]; 22707e43b41SHong Zhang __syncthreads(); 22807e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 22907e43b41SHong Zhang __syncthreads(); 23007e43b41SHong Zhang if (threadIdx.y < 1) { 23107e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 23207e43b41SHong Zhang y[row] = shared[threadIdx.x]; 23307e43b41SHong Zhang } 23407e43b41SHong Zhang } 23507e43b41SHong Zhang } 23607e43b41SHong Zhang 23707e43b41SHong Zhang __global__ void matmult_seqsell_tiled_kernel2(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, PetscScalar *y) 23807e43b41SHong Zhang { 23907e43b41SHong Zhang __shared__ MatScalar shared[512]; 24007e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 24107e43b41SHong Zhang /* multiple threads per row. */ 24207e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 24307e43b41SHong Zhang if (row < nrows) { 24407e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 24507e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 24607e43b41SHong Zhang 24707e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 24807e43b41SHong 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]]; 24907e43b41SHong Zhang __syncthreads(); 25007e43b41SHong Zhang if (threadIdx.y < 1) { 25107e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 25207e43b41SHong Zhang y[row] = shared[threadIdx.x]; 25307e43b41SHong Zhang } 25407e43b41SHong Zhang } 25507e43b41SHong Zhang } 25607e43b41SHong Zhang 25707e43b41SHong Zhang __global__ void matmultadd_seqsell_tiled_kernel6(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z) 25807e43b41SHong Zhang { 25907e43b41SHong Zhang __shared__ MatScalar shared[512]; 26007e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 26107e43b41SHong Zhang /* multiple threads per row. */ 26207e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 26307e43b41SHong Zhang if (row < nrows) { 26407e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 26507e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 26607e43b41SHong Zhang 26707e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 26807e43b41SHong 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]]; 26907e43b41SHong Zhang __syncthreads(); 27007e43b41SHong Zhang if (threadIdx.y < 16) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 16) * blockDim.x + threadIdx.x]; } 27107e43b41SHong Zhang __syncthreads(); 27207e43b41SHong Zhang if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; } 27307e43b41SHong Zhang __syncthreads(); 27407e43b41SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 27507e43b41SHong Zhang __syncthreads(); 27607e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 27707e43b41SHong Zhang __syncthreads(); 27807e43b41SHong Zhang if (threadIdx.y < 1) { 27907e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 2802d1451d4SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 2812d1451d4SHong Zhang } 2822d1451d4SHong Zhang } 2832d1451d4SHong Zhang } 28407e43b41SHong Zhang 28507e43b41SHong Zhang __global__ void matmultadd_seqsell_tiled_kernel5(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z) 28607e43b41SHong Zhang { 28707e43b41SHong Zhang __shared__ MatScalar shared[512]; 28807e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 28907e43b41SHong Zhang /* multiple threads per row. */ 29007e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 29107e43b41SHong Zhang if (row < nrows) { 29207e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 29307e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 29407e43b41SHong Zhang 29507e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 29607e43b41SHong 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]]; 29707e43b41SHong Zhang __syncthreads(); 29807e43b41SHong Zhang if (threadIdx.y < 8) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 8) * blockDim.x + threadIdx.x]; } 29907e43b41SHong Zhang __syncthreads(); 30007e43b41SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 30107e43b41SHong Zhang __syncthreads(); 30207e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 30307e43b41SHong Zhang __syncthreads(); 30407e43b41SHong Zhang if (threadIdx.y < 1) { 30507e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 30607e43b41SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 30707e43b41SHong Zhang } 30807e43b41SHong Zhang } 30907e43b41SHong Zhang } 31007e43b41SHong Zhang 31107e43b41SHong Zhang __global__ void matmultadd_seqsell_tiled_kernel4(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z) 31207e43b41SHong Zhang { 31307e43b41SHong Zhang __shared__ MatScalar shared[512]; 31407e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 31507e43b41SHong Zhang /* multiple threads per row. */ 31607e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 31707e43b41SHong Zhang if (row < nrows) { 31807e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 31907e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 32007e43b41SHong Zhang 32107e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 32207e43b41SHong 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]]; 32307e43b41SHong Zhang __syncthreads(); 32407e43b41SHong Zhang if (threadIdx.y < 4) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 4) * blockDim.x + threadIdx.x]; } 32507e43b41SHong Zhang __syncthreads(); 32607e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 32707e43b41SHong Zhang __syncthreads(); 32807e43b41SHong Zhang if (threadIdx.y < 1) { 32907e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 33007e43b41SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 33107e43b41SHong Zhang } 33207e43b41SHong Zhang } 33307e43b41SHong Zhang } 33407e43b41SHong Zhang 33507e43b41SHong Zhang __global__ void matmultadd_seqsell_tiled_kernel3(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z) 33607e43b41SHong Zhang { 33707e43b41SHong Zhang __shared__ MatScalar shared[512]; 33807e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 33907e43b41SHong Zhang /* multiple threads per row. */ 34007e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 34107e43b41SHong Zhang if (row < nrows) { 34207e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 34307e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 34407e43b41SHong Zhang 34507e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 34607e43b41SHong 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]]; 34707e43b41SHong Zhang __syncthreads(); 34807e43b41SHong Zhang if (threadIdx.y < 2) { shared[threadIdx.y * blockDim.x + threadIdx.x] += shared[(threadIdx.y + 2) * blockDim.x + threadIdx.x]; } 34907e43b41SHong Zhang __syncthreads(); 35007e43b41SHong Zhang if (threadIdx.y < 1) { 35107e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 35207e43b41SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 35307e43b41SHong Zhang } 35407e43b41SHong Zhang } 35507e43b41SHong Zhang } 35607e43b41SHong Zhang 35707e43b41SHong Zhang __global__ void matmultadd_seqsell_tiled_kernel2(PetscInt nrows, PetscInt totalslices, const PetscInt *acolidx, const MatScalar *aval, const PetscInt *sliidx, const PetscScalar *x, const PetscScalar *y, PetscScalar *z) 35807e43b41SHong Zhang { 35907e43b41SHong Zhang __shared__ MatScalar shared[512]; 36007e43b41SHong Zhang PetscInt i, row, slice_id, row_in_slice; 36107e43b41SHong Zhang /* multiple threads per row. */ 36207e43b41SHong Zhang row = blockIdx.x * blockDim.x + threadIdx.x; 36307e43b41SHong Zhang if (row < nrows) { 36407e43b41SHong Zhang slice_id = row / SLICE_HEIGHT; 36507e43b41SHong Zhang row_in_slice = row % SLICE_HEIGHT; 36607e43b41SHong Zhang 36707e43b41SHong Zhang shared[threadIdx.y * blockDim.x + threadIdx.x] = 0.0; 36807e43b41SHong 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]]; 36907e43b41SHong Zhang __syncthreads(); 37007e43b41SHong Zhang if (threadIdx.y < 1) { 37107e43b41SHong Zhang shared[threadIdx.x] += shared[blockDim.x + threadIdx.x]; 37207e43b41SHong Zhang z[row] = y[row] + shared[threadIdx.x]; 37307e43b41SHong Zhang } 37407e43b41SHong Zhang } 3752d1451d4SHong Zhang } 3762d1451d4SHong Zhang 3772d1451d4SHong Zhang PetscErrorCode MatMult_SeqSELLCUDA(Mat A, Vec xx, Vec yy) 3782d1451d4SHong Zhang { 3792d1451d4SHong Zhang Mat_SeqSELL *a = (Mat_SeqSELL *)A->data; 3802d1451d4SHong Zhang Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr; 3812d1451d4SHong Zhang PetscScalar *y; 3822d1451d4SHong Zhang const PetscScalar *x; 383*4e58db63SHong Zhang PetscInt totalslices = a->totalslices, nrows = A->rmap->n, sliceheight = a->sliceheight; 3842d1451d4SHong Zhang MatScalar *aval; 3852d1451d4SHong Zhang PetscInt *acolidx; 3862d1451d4SHong Zhang PetscInt *sliidx; 38707e43b41SHong Zhang PetscInt nblocks, blocksize = 512; /* blocksize must be multiple of SLICE_HEIGHT*32 */ 38807e43b41SHong Zhang dim3 block2(256, 2), block4(128, 4), block8(64, 8), block16(32, 16), block32(16, 32); 3892d1451d4SHong Zhang 3902d1451d4SHong Zhang PetscFunctionBegin; 391*4e58db63SHong 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); 3922d1451d4SHong Zhang PetscCall(MatSeqSELLCUDACopyToGPU(A)); 3932d1451d4SHong Zhang /* cudastruct may not be available until MatSeqSELLCUDACopyToGPU() is called */ 3942d1451d4SHong Zhang aval = cudastruct->val; 3952d1451d4SHong Zhang acolidx = cudastruct->colidx; 3962d1451d4SHong Zhang sliidx = cudastruct->sliidx; 3972d1451d4SHong Zhang 3982d1451d4SHong Zhang PetscCall(VecCUDAGetArrayRead(xx, &x)); 3992d1451d4SHong Zhang PetscCall(VecCUDAGetArrayWrite(yy, &y)); 4002d1451d4SHong Zhang PetscCall(PetscLogGpuTimeBegin()); 40107e43b41SHong Zhang 40207e43b41SHong Zhang switch (cudastruct->kernelchoice) { 40307e43b41SHong Zhang case 9: 404*4e58db63SHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 405*4e58db63SHong Zhang if (cudastruct->blocky == 2) { 406*4e58db63SHong Zhang matmult_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 407*4e58db63SHong Zhang } else if (cudastruct->blocky == 4) { 408*4e58db63SHong Zhang matmult_seqsell_tiled_kernel9<4><<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 409*4e58db63SHong Zhang } else if (cudastruct->blocky == 8) { 410*4e58db63SHong Zhang matmult_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 411*4e58db63SHong Zhang } else if (cudastruct->blocky == 16) { 412*4e58db63SHong Zhang matmult_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 413*4e58db63SHong Zhang } else if (cudastruct->blocky == 32) { 414*4e58db63SHong Zhang matmult_seqsell_tiled_kernel9<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 41507e43b41SHong Zhang } 41607e43b41SHong Zhang break; 41707e43b41SHong Zhang case 7: 418*4e58db63SHong Zhang nblocks = 1 + (nrows - 1) / (2 * sliceheight); 419*4e58db63SHong Zhang if (cudastruct->blocky == 2) { 420*4e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 421*4e58db63SHong Zhang } else if (cudastruct->blocky == 4) { 422*4e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 4)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 423*4e58db63SHong Zhang } else if (cudastruct->blocky == 8) { 424*4e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 425*4e58db63SHong Zhang } else if (cudastruct->blocky == 16) { 426*4e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 427*4e58db63SHong Zhang } else if (cudastruct->blocky == 32) { 428*4e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 429*4e58db63SHong Zhang } else { 430*4e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 431*4e58db63SHong Zhang } 43207e43b41SHong Zhang break; 43307e43b41SHong Zhang case 6: 43407e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 32); /* 1 slice per block if blocksize=512 */ 43507e43b41SHong Zhang matmult_seqsell_tiled_kernel6<<<nblocks, block32>>>(nrows, totalslices, acolidx, aval, sliidx, x, y); 43607e43b41SHong Zhang break; 43707e43b41SHong Zhang case 5: 43807e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 16); /* 2 slices per block if blocksize=512*/ 43907e43b41SHong Zhang matmult_seqsell_tiled_kernel5<<<nblocks, block16>>>(nrows, totalslices, acolidx, aval, sliidx, x, y); 44007e43b41SHong Zhang break; 44107e43b41SHong Zhang case 4: 44207e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 8); /* 4 slices per block if blocksize=512 */ 44307e43b41SHong Zhang matmult_seqsell_tiled_kernel4<<<nblocks, block8>>>(nrows, totalslices, acolidx, aval, sliidx, x, y); 44407e43b41SHong Zhang break; 44507e43b41SHong Zhang case 3: 44607e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 4); /* 8 slices per block if blocksize=512 */ 44707e43b41SHong Zhang matmult_seqsell_tiled_kernel3<<<nblocks, block4>>>(nrows, totalslices, acolidx, aval, sliidx, x, y); 44807e43b41SHong Zhang break; 44907e43b41SHong Zhang case 2: /* 16 slices per block if blocksize=512 */ 45007e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 2); 45107e43b41SHong Zhang matmult_seqsell_tiled_kernel2<<<nblocks, block2>>>(nrows, totalslices, acolidx, aval, sliidx, x, y); 45207e43b41SHong Zhang break; 45307e43b41SHong Zhang case 1: /* 32 slices per block if blocksize=512 */ 45407e43b41SHong Zhang nblocks = 1 + (nrows - 1) / blocksize; 455*4e58db63SHong Zhang matmult_seqsell_basic_kernel<<<nblocks, blocksize>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 45607e43b41SHong Zhang break; 45707e43b41SHong Zhang case 0: 458*4e58db63SHong Zhang if (sliceheight * a->maxslicewidth > 20800) { /* important threshold */ 459*4e58db63SHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 460*4e58db63SHong Zhang matmult_seqsell_tiled_kernel9<32><<<nblocks, dim3(32, 32)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 4612d1451d4SHong Zhang } else { 462*4e58db63SHong Zhang PetscInt avgslicesize = sliceheight * a->avgslicewidth; 463*4e58db63SHong Zhang if (avgslicesize <= 96) { 464*4e58db63SHong Zhang nblocks = 1 + (nrows - 1) / (2 * sliceheight); /* two slices per block */ 465*4e58db63SHong Zhang matmult_seqsell_tiled_kernel7<<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 466*4e58db63SHong Zhang } else if (avgslicesize <= 432) { 467*4e58db63SHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 468*4e58db63SHong Zhang matmult_seqsell_tiled_kernel9<2><<<nblocks, dim3(32, 2)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 469*4e58db63SHong Zhang } else if (avgslicesize <= 2400) { 470*4e58db63SHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 471*4e58db63SHong Zhang matmult_seqsell_tiled_kernel9<8><<<nblocks, dim3(32, 8)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 4722d1451d4SHong Zhang } else { 473*4e58db63SHong Zhang nblocks = 1 + (nrows - 1) / sliceheight; 474*4e58db63SHong Zhang matmult_seqsell_tiled_kernel9<16><<<nblocks, dim3(32, 16)>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y); 4752d1451d4SHong Zhang } 4762d1451d4SHong Zhang } 47707e43b41SHong Zhang break; 47807e43b41SHong Zhang } 4792d1451d4SHong Zhang PetscCall(PetscLogGpuTimeEnd()); 4802d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayRead(xx, &x)); 4812d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayWrite(yy, &y)); 4822d1451d4SHong Zhang PetscCall(PetscLogGpuFlops(2.0 * a->nz - a->nonzerorowcnt)); 4832d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 4842d1451d4SHong Zhang } 4852d1451d4SHong Zhang 4862d1451d4SHong Zhang PetscErrorCode MatMultAdd_SeqSELLCUDA(Mat A, Vec xx, Vec yy, Vec zz) 4872d1451d4SHong Zhang { 4882d1451d4SHong Zhang Mat_SeqSELL *a = (Mat_SeqSELL *)A->data; 4892d1451d4SHong Zhang Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr; 4902d1451d4SHong Zhang PetscScalar *z; 4912d1451d4SHong Zhang const PetscScalar *y, *x; 492*4e58db63SHong Zhang PetscInt totalslices = a->totalslices, nrows = A->rmap->n, sliceheight = a->sliceheight; 4932d1451d4SHong Zhang MatScalar *aval = cudastruct->val; 4942d1451d4SHong Zhang PetscInt *acolidx = cudastruct->colidx; 4952d1451d4SHong Zhang PetscInt *sliidx = cudastruct->sliidx; 4962d1451d4SHong Zhang 4972d1451d4SHong Zhang PetscFunctionBegin; 498*4e58db63SHong Zhang PetscCheck(sliceheight == 16, PETSC_COMM_SELF, PETSC_ERR_SUP, "The kernel requires a slice height of 16, but the input matrix has a slice height of %" PetscInt_FMT, sliceheight); 4992d1451d4SHong Zhang PetscCall(MatSeqSELLCUDACopyToGPU(A)); 5002d1451d4SHong Zhang if (a->nz) { 50107e43b41SHong Zhang PetscInt nblocks, blocksize = 512; 50207e43b41SHong Zhang dim3 block2(256, 2), block4(128, 4), block8(64, 8), block16(32, 16), block32(16, 32); 5032d1451d4SHong Zhang PetscCall(VecCUDAGetArrayRead(xx, &x)); 5042d1451d4SHong Zhang PetscCall(VecCUDAGetArrayRead(yy, &y)); 5052d1451d4SHong Zhang PetscCall(VecCUDAGetArrayWrite(zz, &z)); 5062d1451d4SHong Zhang PetscCall(PetscLogGpuTimeBegin()); 50707e43b41SHong Zhang 50807e43b41SHong Zhang switch (cudastruct->kernelchoice) { 50907e43b41SHong Zhang case 6: 51007e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 32); 51107e43b41SHong Zhang matmultadd_seqsell_tiled_kernel6<<<nblocks, block32>>>(nrows, totalslices, acolidx, aval, sliidx, x, y, z); 51207e43b41SHong Zhang break; 51307e43b41SHong Zhang case 5: 51407e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 16); 51507e43b41SHong Zhang matmultadd_seqsell_tiled_kernel5<<<nblocks, block16>>>(nrows, totalslices, acolidx, aval, sliidx, x, y, z); 51607e43b41SHong Zhang break; 51707e43b41SHong Zhang case 4: 51807e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 8); 51907e43b41SHong Zhang matmultadd_seqsell_tiled_kernel4<<<nblocks, block8>>>(nrows, totalslices, acolidx, aval, sliidx, x, y, z); 52007e43b41SHong Zhang break; 52107e43b41SHong Zhang case 3: 52207e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 4); 52307e43b41SHong Zhang matmultadd_seqsell_tiled_kernel3<<<nblocks, block4>>>(nrows, totalslices, acolidx, aval, sliidx, x, y, z); 52407e43b41SHong Zhang break; 52507e43b41SHong Zhang case 2: 52607e43b41SHong Zhang nblocks = 1 + (nrows - 1) / (blocksize / 2); 52707e43b41SHong Zhang matmultadd_seqsell_tiled_kernel2<<<nblocks, block2>>>(nrows, totalslices, acolidx, aval, sliidx, x, y, z); 52807e43b41SHong Zhang break; 52907e43b41SHong Zhang case 1: 53007e43b41SHong Zhang nblocks = 1 + (nrows - 1) / blocksize; 531*4e58db63SHong Zhang matmultadd_seqsell_basic_kernel<<<nblocks, blocksize>>>(nrows, sliceheight, acolidx, aval, sliidx, x, y, z); 53207e43b41SHong Zhang break; 53307e43b41SHong Zhang case 0: /* TODO */ 53407e43b41SHong Zhang break; 5352d1451d4SHong Zhang } 5362d1451d4SHong Zhang PetscCall(PetscLogGpuTimeEnd()); 5372d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayRead(xx, &x)); 5382d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayRead(yy, &y)); 5392d1451d4SHong Zhang PetscCall(VecCUDARestoreArrayWrite(zz, &z)); 5402d1451d4SHong Zhang PetscCall(PetscLogGpuFlops(2.0 * a->nz)); 5412d1451d4SHong Zhang } else { 5422d1451d4SHong Zhang PetscCall(VecCopy(yy, zz)); 5432d1451d4SHong Zhang } 5442d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 5452d1451d4SHong Zhang } 5462d1451d4SHong Zhang 5472d1451d4SHong Zhang static PetscErrorCode MatSetFromOptions_SeqSELLCUDA(Mat A, PetscOptionItems *PetscOptionsObject) 5482d1451d4SHong Zhang { 54907e43b41SHong Zhang Mat_SeqSELLCUDA *cudastruct = (Mat_SeqSELLCUDA *)A->spptr; 550*4e58db63SHong Zhang PetscInt kernel, blocky; 55107e43b41SHong Zhang PetscBool flg; 55207e43b41SHong Zhang 5532d1451d4SHong Zhang PetscFunctionBegin; 5542d1451d4SHong Zhang PetscOptionsHeadBegin(PetscOptionsObject, "SeqSELLCUDA options"); 55507e43b41SHong Zhang PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_sell_spmv_cuda_kernel", &kernel, &flg)); 55607e43b41SHong Zhang if (flg) { 55707e43b41SHong 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); 55807e43b41SHong Zhang cudastruct->kernelchoice = kernel; 55907e43b41SHong Zhang } 560*4e58db63SHong Zhang PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_sell_spmv_cuda_blocky", &blocky, &flg)); 561*4e58db63SHong Zhang if (flg) { 562*4e58db63SHong 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}", kernel); 563*4e58db63SHong Zhang cudastruct->blocky = blocky; 564*4e58db63SHong Zhang } 5652d1451d4SHong Zhang PetscOptionsHeadEnd(); 5662d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 5672d1451d4SHong Zhang } 5682d1451d4SHong Zhang 56907e43b41SHong Zhang PETSC_INTERN PetscErrorCode MatAssemblyEnd_SpMV_Preprocessing_Private(Mat A) 57007e43b41SHong Zhang { 57107e43b41SHong Zhang Mat_SeqSELL *a = (Mat_SeqSELL *)A->data; 57207e43b41SHong Zhang 57307e43b41SHong Zhang PetscCall(MatSeqSELLGetAvgSliceWidth(A, &a->avgslicewidth)); 57407e43b41SHong Zhang PetscCall(MatSeqSELLGetMaxSliceWidth(A, &a->maxslicewidth)); 57507e43b41SHong Zhang PetscCall(MatSeqSELLGetFillRatio(A, &a->fillratio)); 57607e43b41SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 57707e43b41SHong Zhang } 57807e43b41SHong Zhang 5792d1451d4SHong Zhang static PetscErrorCode MatAssemblyEnd_SeqSELLCUDA(Mat A, MatAssemblyType mode) 5802d1451d4SHong Zhang { 5812d1451d4SHong Zhang PetscFunctionBegin; 5822d1451d4SHong Zhang PetscCall(MatAssemblyEnd_SeqSELL(A, mode)); 58307e43b41SHong Zhang PetscCall(MatAssemblyEnd_SpMV_Preprocessing_Private(A)); 5842d1451d4SHong Zhang if (mode == MAT_FLUSH_ASSEMBLY) PetscFunctionReturn(PETSC_SUCCESS); 5852d1451d4SHong Zhang if (A->factortype == MAT_FACTOR_NONE) { PetscCall(MatSeqSELLCUDACopyToGPU(A)); } 5862d1451d4SHong Zhang A->ops->mult = MatMult_SeqSELLCUDA; 5872d1451d4SHong Zhang A->ops->multadd = MatMultAdd_SeqSELLCUDA; 5882d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 5892d1451d4SHong Zhang } 5902d1451d4SHong Zhang 5912d1451d4SHong Zhang static PetscErrorCode MatDestroy_SeqSELLCUDA(Mat A) 5922d1451d4SHong Zhang { 5932d1451d4SHong Zhang PetscFunctionBegin; 5942d1451d4SHong Zhang if (A->factortype == MAT_FACTOR_NONE) { 5952d1451d4SHong Zhang if (A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) { PetscCall(MatSeqSELLCUDA_Destroy((Mat_SeqSELLCUDA **)&A->spptr)); } 5962d1451d4SHong Zhang } 5972d1451d4SHong Zhang PetscCall(MatDestroy_SeqSELL(A)); 5982d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 5992d1451d4SHong Zhang } 6002d1451d4SHong Zhang 6012d1451d4SHong Zhang PETSC_INTERN PetscErrorCode MatConvert_SeqSELL_SeqSELLCUDA(Mat); 6022d1451d4SHong Zhang static PetscErrorCode MatDuplicate_SeqSELLCUDA(Mat A, MatDuplicateOption cpvalues, Mat *B) 6032d1451d4SHong Zhang { 6042d1451d4SHong Zhang PetscFunctionBegin; 6052d1451d4SHong Zhang PetscCall(MatDuplicate_SeqSELL(A, cpvalues, B)); 6062d1451d4SHong Zhang PetscCall(MatConvert_SeqSELL_SeqSELLCUDA(*B)); 6072d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 6082d1451d4SHong Zhang } 6092d1451d4SHong Zhang 6102d1451d4SHong Zhang PETSC_EXTERN PetscErrorCode MatConvert_SeqSELL_SeqSELLCUDA(Mat B) 6112d1451d4SHong Zhang { 6122d1451d4SHong Zhang Mat_SeqSELLCUDA *cudastruct; 6132d1451d4SHong Zhang 6142d1451d4SHong Zhang PetscFunctionBegin; 6152d1451d4SHong Zhang PetscCall(PetscFree(B->defaultvectype)); 6162d1451d4SHong Zhang PetscCall(PetscStrallocpy(VECCUDA, &B->defaultvectype)); 6172d1451d4SHong Zhang 6182d1451d4SHong Zhang if (!B->spptr) { 6192d1451d4SHong Zhang if (B->factortype == MAT_FACTOR_NONE) { 6202d1451d4SHong Zhang PetscCall(PetscNew(&cudastruct)); 6212d1451d4SHong Zhang B->spptr = cudastruct; 6222d1451d4SHong Zhang } 6232d1451d4SHong Zhang } 6242d1451d4SHong Zhang 6252d1451d4SHong Zhang B->ops->assemblyend = MatAssemblyEnd_SeqSELLCUDA; 6262d1451d4SHong Zhang B->ops->destroy = MatDestroy_SeqSELLCUDA; 6272d1451d4SHong Zhang B->ops->setfromoptions = MatSetFromOptions_SeqSELLCUDA; 6282d1451d4SHong Zhang B->ops->mult = MatMult_SeqSELLCUDA; 6292d1451d4SHong Zhang B->ops->multadd = MatMultAdd_SeqSELLCUDA; 6302d1451d4SHong Zhang B->ops->duplicate = MatDuplicate_SeqSELLCUDA; 6312d1451d4SHong Zhang 63207e43b41SHong Zhang /* No need to assemble SeqSELL, but need to do the preprocessing for SpMV */ 63307e43b41SHong Zhang PetscCall(MatAssemblyEnd_SpMV_Preprocessing_Private(B)); 63407e43b41SHong Zhang 6352d1451d4SHong Zhang PetscCall(PetscObjectChangeTypeName((PetscObject)B, MATSEQSELLCUDA)); 6362d1451d4SHong Zhang B->offloadmask = PETSC_OFFLOAD_UNALLOCATED; 6372d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 6382d1451d4SHong Zhang } 6392d1451d4SHong Zhang 6402d1451d4SHong Zhang PETSC_EXTERN PetscErrorCode MatCreate_SeqSELLCUDA(Mat B) 6412d1451d4SHong Zhang { 6422d1451d4SHong Zhang PetscFunctionBegin; 6432d1451d4SHong Zhang PetscCall(MatCreate_SeqSELL(B)); 6442d1451d4SHong Zhang PetscCall(MatConvert_SeqSELL_SeqSELLCUDA(B)); 6452d1451d4SHong Zhang PetscFunctionReturn(PETSC_SUCCESS); 6462d1451d4SHong Zhang } 647