xref: /petsc/src/mat/impls/hypre/cuda/hypre1.cu (revision 3573c4fa955bc17c426001724e09bd081ebd8244)
1*78923ab9SStefano Zampini #include <petsc/private/petschypre.h>
2a32e9c99SJunchao Zhang #include <petscdevice_cuda.h>
3a32e9c99SJunchao Zhang #include <../src/mat/impls/hypre/mhypre_kernels.hpp>
4a32e9c99SJunchao Zhang 
MatZeroRows_CUDA(PetscInt n,const PetscInt rows[],const HYPRE_Int i[],const HYPRE_Int j[],HYPRE_Complex a[],HYPRE_Complex diag)5a32e9c99SJunchao Zhang PetscErrorCode MatZeroRows_CUDA(PetscInt n, const PetscInt rows[], const HYPRE_Int i[], const HYPRE_Int j[], HYPRE_Complex a[], HYPRE_Complex diag)
6a32e9c99SJunchao Zhang {
7a32e9c99SJunchao Zhang   const PetscInt blkDimX = 16, blkDimY = 32;
8a32e9c99SJunchao Zhang   PetscInt       gridDimX = (n + blkDimX - 1) / blkDimX;
9a32e9c99SJunchao Zhang   cudaStream_t   stream;
10a32e9c99SJunchao Zhang 
11a32e9c99SJunchao Zhang   PetscFunctionBegin;
12a32e9c99SJunchao Zhang   if (!n) PetscFunctionReturn(PETSC_SUCCESS);
13a32e9c99SJunchao Zhang   PetscCall(PetscGetCurrentCUDAStream(&stream));
14a32e9c99SJunchao Zhang   ZeroRows<<<dim3(gridDimX, 1), dim3(blkDimX, blkDimY), 0, stream>>>(n, rows, i, j, a, diag);
15a32e9c99SJunchao Zhang   PetscCallCUDA(cudaGetLastError());
16a32e9c99SJunchao Zhang   PetscFunctionReturn(PETSC_SUCCESS);
17a32e9c99SJunchao Zhang }
18