1a4963045SJacob Faibussowitsch #pragma once 292896123SJunchao Zhang #include <petsc_kokkos.hpp> 3c0c276a7Ssdargavi #include <petscmat_kokkos.hpp> 4e907feaaSJunchao Zhang #include <petsc/private/kokkosimpl.hpp> 5394ed5ebSJunchao Zhang #include <../src/mat/impls/aij/seq/aij.h> 642550becSJunchao Zhang #include <KokkosSparse_CrsMatrix.hpp> 742550becSJunchao Zhang #include <KokkosSparse_spiluk.hpp> 87b8d4ba6SJunchao Zhang #include <string> 97b8d4ba6SJunchao Zhang 107b8d4ba6SJunchao Zhang namespace 117b8d4ba6SJunchao Zhang { 127b8d4ba6SJunchao Zhang PETSC_NODISCARD inline decltype(auto) NoInit(std::string label) 137b8d4ba6SJunchao Zhang { 147b8d4ba6SJunchao Zhang return Kokkos::view_alloc(Kokkos::WithoutInitializing, std::move(label)); 157b8d4ba6SJunchao Zhang } 167b8d4ba6SJunchao Zhang } // namespace 1742550becSJunchao Zhang 1842550becSJunchao Zhang using MatRowMapType = PetscInt; 1942550becSJunchao Zhang using MatColIdxType = PetscInt; 2042550becSJunchao Zhang using MatScalarType = PetscScalar; 2142550becSJunchao Zhang 229371c9d4SSatish Balay template <class MemorySpace> 239371c9d4SSatish Balay using KokkosCsrMatrixType = typename KokkosSparse::CrsMatrix<MatScalarType, MatColIdxType, MemorySpace, void /* MemoryTraits */, MatRowMapType>; 249371c9d4SSatish Balay template <class MemorySpace> 259371c9d4SSatish Balay using KokkosCsrGraphType = typename KokkosCsrMatrixType<MemorySpace>::staticcrsgraph_type; 2642550becSJunchao Zhang 2742550becSJunchao Zhang using KokkosCsrGraph = KokkosCsrGraphType<DefaultMemorySpace>; 2845402d8aSJunchao Zhang using KokkosCsrGraphHost = KokkosCsrGraphType<HostMirrorMemorySpace>; 2942550becSJunchao Zhang 3042550becSJunchao Zhang using KokkosCsrMatrix = KokkosCsrMatrixType<DefaultMemorySpace>; 3145402d8aSJunchao Zhang using KokkosCsrMatrixHost = KokkosCsrMatrixType<HostMirrorMemorySpace>; 3242550becSJunchao Zhang 3342550becSJunchao Zhang using MatRowMapKokkosView = KokkosCsrGraph::row_map_type::non_const_type; 3442550becSJunchao Zhang using MatColIdxKokkosView = KokkosCsrGraph::entries_type::non_const_type; 3542550becSJunchao Zhang using MatScalarKokkosView = KokkosCsrMatrix::values_type::non_const_type; 3642550becSJunchao Zhang 3742550becSJunchao Zhang using MatRowMapKokkosViewHost = KokkosCsrGraphHost::row_map_type::non_const_type; 3842550becSJunchao Zhang using MatColIdxKokkosViewHost = KokkosCsrGraphHost::entries_type::non_const_type; 3942550becSJunchao Zhang using MatScalarKokkosViewHost = KokkosCsrMatrixHost::values_type::non_const_type; 4042550becSJunchao Zhang 4142550becSJunchao Zhang using ConstMatRowMapKokkosView = KokkosCsrGraph::row_map_type::const_type; 4242550becSJunchao Zhang using ConstMatColIdxKokkosView = KokkosCsrGraph::entries_type::const_type; 4342550becSJunchao Zhang using ConstMatScalarKokkosView = KokkosCsrMatrix::values_type::const_type; 4442550becSJunchao Zhang 4542550becSJunchao Zhang using ConstMatRowMapKokkosViewHost = KokkosCsrGraphHost::row_map_type::const_type; 4642550becSJunchao Zhang using ConstMatColIdxKokkosViewHost = KokkosCsrGraphHost::entries_type::const_type; 4742550becSJunchao Zhang using ConstMatScalarKokkosViewHost = KokkosCsrMatrixHost::values_type::const_type; 4842550becSJunchao Zhang 4942550becSJunchao Zhang using MatRowMapKokkosDualView = Kokkos::DualView<MatRowMapType *>; 5042550becSJunchao Zhang using MatColIdxKokkosDualView = Kokkos::DualView<MatColIdxType *>; 5142550becSJunchao Zhang using MatScalarKokkosDualView = Kokkos::DualView<MatScalarType *>; 5242550becSJunchao Zhang 5342550becSJunchao Zhang using KernelHandle = KokkosKernels::Experimental::KokkosKernelsHandle<MatRowMapType, MatColIdxType, MatScalarType, DefaultExecutionSpace, DefaultMemorySpace, DefaultMemorySpace>; 5442550becSJunchao Zhang 5542550becSJunchao Zhang using KokkosTeamMemberType = Kokkos::TeamPolicy<DefaultExecutionSpace>::member_type; 5642550becSJunchao Zhang 5742550becSJunchao Zhang /* For mat->spptr of a factorized matrix */ 5842550becSJunchao Zhang struct Mat_SeqAIJKokkosTriFactors { 5942550becSJunchao Zhang MatRowMapKokkosView iL_d, iU_d, iLt_d, iUt_d; /* rowmap for L, U, L^t, U^t of A=LU */ 6042550becSJunchao Zhang MatColIdxKokkosView jL_d, jU_d, jLt_d, jUt_d; /* column ids */ 6142550becSJunchao Zhang MatScalarKokkosView aL_d, aU_d, aLt_d, aUt_d; /* matrix values */ 62aac854edSJunchao Zhang KernelHandle kh, khL, khU, khLt, khUt; /* Kernel handles for ILU factorization of A, and TRSV of L, U, L^t, U^t */ 63aac854edSJunchao Zhang PetscScalarKokkosView workVector; 6442550becSJunchao Zhang PetscBool transpose_updated; /* Are L^T, U^T updated wrt L, U*/ 6542550becSJunchao Zhang PetscBool sptrsv_symbolic_completed; /* Have we completed the symbolic solve for L and U */ 6642550becSJunchao Zhang 67aac854edSJunchao Zhang MatRowMapKokkosViewHost iL_h, iU_h, iLt_h, iUt_h; // temp. buffers when we do factorization with PETSc on host. We copy L, U to these buffers and then copy to device 68aac854edSJunchao Zhang MatColIdxKokkosViewHost jL_h, jU_h, jLt_h, jUt_h; 69aac854edSJunchao Zhang MatScalarKokkosViewHost aL_h, aU_h, aLt_h, aUt_h, D_h; // D is for LDLT factorization 70aac854edSJunchao Zhang MatScalarKokkosView D_d; 71aac854edSJunchao Zhang Mat L, U, Lt, Ut; // MATSEQAIJ on host if needed. Their arrays are alias to (iL_h, jL_h, aL_h), (iU_h, jU_h, aU_h) and their transpose. 72aac854edSJunchao Zhang // MatTranspose() on host might be faster than KK's csr transpose on device. 73aac854edSJunchao Zhang 74aac854edSJunchao Zhang PetscIntKokkosView rowperm, colperm; // row permutation and column permutation 75aac854edSJunchao Zhang 76aac854edSJunchao Zhang Mat_SeqAIJKokkosTriFactors(PetscInt n) : workVector("workVector", n) 77aac854edSJunchao Zhang { 78aac854edSJunchao Zhang L = U = Lt = Ut = nullptr; 79aac854edSJunchao Zhang transpose_updated = sptrsv_symbolic_completed = PETSC_FALSE; 80aac854edSJunchao Zhang } 8142550becSJunchao Zhang 8242550becSJunchao Zhang ~Mat_SeqAIJKokkosTriFactors() { Destroy(); } 8342550becSJunchao Zhang 84d71ae5a4SJacob Faibussowitsch void Destroy() 85d71ae5a4SJacob Faibussowitsch { 86aac854edSJunchao Zhang PetscFunctionBeginUser; 8742550becSJunchao Zhang kh.destroy_spiluk_handle(); 8842550becSJunchao Zhang khL.destroy_sptrsv_handle(); 8942550becSJunchao Zhang khU.destroy_sptrsv_handle(); 9042550becSJunchao Zhang khLt.destroy_sptrsv_handle(); 9142550becSJunchao Zhang khUt.destroy_sptrsv_handle(); 92aac854edSJunchao Zhang PetscCallVoid(MatDestroy(&L)); 93aac854edSJunchao Zhang PetscCallVoid(MatDestroy(&U)); 94aac854edSJunchao Zhang PetscCallVoid(MatDestroy(&Lt)); 95aac854edSJunchao Zhang PetscCallVoid(MatDestroy(&Ut)); 96aac854edSJunchao Zhang L = U = Lt = Ut = nullptr; 9742550becSJunchao Zhang transpose_updated = sptrsv_symbolic_completed = PETSC_FALSE; 98aac854edSJunchao Zhang PetscFunctionReturnVoid(); 9942550becSJunchao Zhang } 10042550becSJunchao Zhang }; 10142550becSJunchao Zhang 10242550becSJunchao Zhang /* For mat->spptr of a regular matrix */ 10342550becSJunchao Zhang struct Mat_SeqAIJKokkos { 10442550becSJunchao Zhang MatRowMapKokkosDualView i_dual; 10542550becSJunchao Zhang MatColIdxKokkosDualView j_dual; 10642550becSJunchao Zhang MatScalarKokkosDualView a_dual; 10742550becSJunchao Zhang 108f78ce678SMark Adams MatRowMapKokkosDualView diag_dual; /* Diagonal pointer, built on demand */ 109f78ce678SMark Adams 11042550becSJunchao Zhang KokkosCsrMatrix csrmat; /* The CSR matrix, used to call KK functions */ 11142550becSJunchao Zhang PetscObjectState nonzerostate; /* State of the nonzero pattern (graph) on device */ 11242550becSJunchao Zhang 11342550becSJunchao Zhang KokkosCsrMatrix csrmatT, csrmatH; /* Transpose and Hermitian of the matrix (built on demand) */ 11442550becSJunchao Zhang PetscBool transpose_updated, hermitian_updated; /* Are At, Ah updated wrt the matrix? */ 1150e3ece09SJunchao Zhang MatRowMapKokkosView transpose_perm; // A permutation array making Ta(i) = Aa(perm(i)), where T = A^t 11642550becSJunchao Zhang 117f4747e26SJunchao Zhang /* Construct a nrows by ncols matrix with given aseq on host. Caller also specifies a nonzero state */ 118f4747e26SJunchao Zhang Mat_SeqAIJKokkos(PetscInt nrows, PetscInt ncols, Mat_SeqAIJ *aseq, PetscObjectState nzstate, PetscBool copyValues = PETSC_TRUE) 119d71ae5a4SJacob Faibussowitsch { 1204df4a32cSJunchao Zhang auto exec = PetscGetKokkosExecutionSpace(); 12192896123SJunchao Zhang 122f4747e26SJunchao Zhang MatScalarKokkosViewHost a_h(aseq->a, aseq->nz); 123f4747e26SJunchao Zhang MatRowMapKokkosViewHost i_h(const_cast<MatRowMapType *>(aseq->i), nrows + 1); 124f4747e26SJunchao Zhang MatColIdxKokkosViewHost j_h(aseq->j, aseq->nz); 125f4747e26SJunchao Zhang MatRowMapKokkosViewHost diag_h(aseq->diag, nrows); 12642550becSJunchao Zhang 12792896123SJunchao Zhang auto a_d = Kokkos::create_mirror_view(Kokkos::WithoutInitializing, exec, a_h); 12892896123SJunchao Zhang auto i_d = Kokkos::create_mirror_view_and_copy(exec, i_h); 12992896123SJunchao Zhang auto j_d = Kokkos::create_mirror_view_and_copy(exec, j_h); 13092896123SJunchao Zhang auto diag_d = Kokkos::create_mirror_view_and_copy(exec, diag_h); 13142550becSJunchao Zhang a_dual = MatScalarKokkosDualView(a_d, a_h); 13242550becSJunchao Zhang i_dual = MatRowMapKokkosDualView(i_d, i_h); 13342550becSJunchao Zhang j_dual = MatColIdxKokkosDualView(j_d, j_h); 134f4747e26SJunchao Zhang diag_dual = MatColIdxKokkosDualView(diag_d, diag_h); 13542550becSJunchao Zhang 13642550becSJunchao Zhang a_dual.modify_host(); /* Since caller provided values on host */ 137*f3d3cd90SZach Atkins if (copyValues) (void)KokkosDualViewSyncDevice(a_dual, exec); 13842550becSJunchao Zhang 13942550becSJunchao Zhang csrmat = KokkosCsrMatrix("csrmat", ncols, a_d, KokkosCsrGraph(j_d, i_d)); 140aac854edSJunchao Zhang Init(nzstate); 14142550becSJunchao Zhang } 14242550becSJunchao Zhang 14342550becSJunchao Zhang /* Construct with a KokkosCsrMatrix. For performance, only i, j are copied to host, but not the matrix values. */ 144d71ae5a4SJacob Faibussowitsch Mat_SeqAIJKokkos(const KokkosCsrMatrix &csr) : csrmat(csr) /* Shallow-copy csr's views to csrmat */ 14542550becSJunchao Zhang { 14642550becSJunchao Zhang auto a_d = csr.values; 14742550becSJunchao Zhang /* Get a non-const version since I don't want to deal with DualView<const T*>, which is not well defined */ 14842550becSJunchao Zhang MatRowMapKokkosView i_d(const_cast<MatRowMapType *>(csr.graph.row_map.data()), csr.graph.row_map.extent(0)); 14942550becSJunchao Zhang auto j_d = csr.graph.entries; 15045402d8aSJunchao Zhang auto a_h = Kokkos::create_mirror_view(Kokkos::WithoutInitializing, HostMirrorMemorySpace(), a_d); 15145402d8aSJunchao Zhang auto i_h = Kokkos::create_mirror_view_and_copy(HostMirrorMemorySpace(), i_d); 15245402d8aSJunchao Zhang auto j_h = Kokkos::create_mirror_view_and_copy(HostMirrorMemorySpace(), j_d); 15342550becSJunchao Zhang 154f4747e26SJunchao Zhang // diag_dual is set until MatAssemblyEnd() where we copy diag from host to device 15542550becSJunchao Zhang a_dual = MatScalarKokkosDualView(a_d, a_h); 15642550becSJunchao Zhang a_dual.modify_device(); /* since we did not copy a_d to a_h, we mark device has the latest data */ 15742550becSJunchao Zhang i_dual = MatRowMapKokkosDualView(i_d, i_h); 15842550becSJunchao Zhang j_dual = MatColIdxKokkosDualView(j_d, j_h); 15942550becSJunchao Zhang Init(); 16042550becSJunchao Zhang } 16142550becSJunchao Zhang 162d71ae5a4SJacob Faibussowitsch Mat_SeqAIJKokkos(PetscInt nrows, PetscInt ncols, PetscInt nnz, MatRowMapKokkosDualView &i, MatColIdxKokkosDualView &j, MatScalarKokkosDualView a) : i_dual(i), j_dual(j), a_dual(a) 163d71ae5a4SJacob Faibussowitsch { 16442550becSJunchao Zhang csrmat = KokkosCsrMatrix("csrmat", nrows, ncols, nnz, a.view_device(), i.view_device(), j.view_device()); 16542550becSJunchao Zhang Init(); 16642550becSJunchao Zhang } 16742550becSJunchao Zhang 16842550becSJunchao Zhang MatScalarType *a_host_data() { return a_dual.view_host().data(); } 16942550becSJunchao Zhang MatRowMapType *i_host_data() { return i_dual.view_host().data(); } 17042550becSJunchao Zhang MatColIdxType *j_host_data() { return j_dual.view_host().data(); } 17142550becSJunchao Zhang 17242550becSJunchao Zhang MatScalarType *a_device_data() { return a_dual.view_device().data(); } 17342550becSJunchao Zhang MatRowMapType *i_device_data() { return i_dual.view_device().data(); } 17442550becSJunchao Zhang MatColIdxType *j_device_data() { return j_dual.view_device().data(); } 17542550becSJunchao Zhang 17642550becSJunchao Zhang MatColIdxType nrows() { return csrmat.numRows(); } 17742550becSJunchao Zhang MatColIdxType ncols() { return csrmat.numCols(); } 17842550becSJunchao Zhang MatRowMapType nnz() { return csrmat.nnz(); } 17942550becSJunchao Zhang 18042550becSJunchao Zhang /* Change the csrmat size to n */ 18142550becSJunchao Zhang void SetColSize(MatColIdxType n) { csrmat = KokkosCsrMatrix("csrmat", n, a_dual.view_device(), csrmat.graph); } 18242550becSJunchao Zhang 183d71ae5a4SJacob Faibussowitsch void SetDiagonal(const MatRowMapType *diag) 184d71ae5a4SJacob Faibussowitsch { 185f78ce678SMark Adams MatRowMapKokkosViewHost diag_h(const_cast<MatRowMapType *>(diag), nrows()); 186f78ce678SMark Adams auto diag_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), diag_h); 187f78ce678SMark Adams diag_dual = MatRowMapKokkosDualView(diag_d, diag_h); 188f78ce678SMark Adams } 189f78ce678SMark Adams 19042550becSJunchao Zhang /* Shared init stuff */ 191aac854edSJunchao Zhang void Init(PetscObjectState nzstate = 0) 192d71ae5a4SJacob Faibussowitsch { 193aac854edSJunchao Zhang nonzerostate = nzstate; 194aac854edSJunchao Zhang transpose_updated = PETSC_FALSE; 195aac854edSJunchao Zhang hermitian_updated = PETSC_FALSE; 19642550becSJunchao Zhang } 19742550becSJunchao Zhang 198d71ae5a4SJacob Faibussowitsch PetscErrorCode DestroyMatTranspose(void) 199d71ae5a4SJacob Faibussowitsch { 20042550becSJunchao Zhang PetscFunctionBegin; 20142550becSJunchao Zhang csrmatT = KokkosCsrMatrix(); /* Overwrite with empty matrices */ 20242550becSJunchao Zhang csrmatH = KokkosCsrMatrix(); 2033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20442550becSJunchao Zhang } 20542550becSJunchao Zhang }; 20642550becSJunchao Zhang 20742550becSJunchao Zhang struct MatProductData_SeqAIJKokkos { 20842550becSJunchao Zhang KernelHandle kh; 20942550becSJunchao Zhang PetscBool reusesym; 21042550becSJunchao Zhang MatProductData_SeqAIJKokkos() : reusesym(PETSC_FALSE) { } 21142550becSJunchao Zhang }; 21242550becSJunchao Zhang 21342550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetSeqAIJKokkosWithCSRMatrix(Mat, Mat_SeqAIJKokkos *); 21442550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatCreateSeqAIJKokkosWithCSRMatrix(MPI_Comm, Mat_SeqAIJKokkos *, Mat *); 21542550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosMergeMats(Mat, Mat, MatReuse, Mat *); 21642550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosSyncDevice(Mat); 2170e3ece09SJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosGetKokkosCsrMatrix(Mat, KokkosCsrMatrix *); 2180e3ece09SJunchao Zhang PETSC_INTERN PetscErrorCode MatCreateSeqAIJKokkosWithKokkosCsrMatrix(MPI_Comm, KokkosCsrMatrix, Mat *); 21942550becSJunchao Zhang PETSC_INTERN PetscErrorCode PrintCsrMatrix(const KokkosCsrMatrix &csrmat); 22042550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJKokkos(Mat, MatType, MatReuse, Mat *); 22142550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosModifyDevice(Mat); 2220e3ece09SJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosGenerateTranspose_Private(Mat, KokkosCsrMatrix *); 223394ed5ebSJunchao Zhang 224394ed5ebSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJGetKokkosView(Mat, MatScalarKokkosView *); 225394ed5ebSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJRestoreKokkosView(Mat, MatScalarKokkosView *); 226394ed5ebSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJGetKokkosViewWrite(Mat, MatScalarKokkosView *); 227394ed5ebSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJRestoreKokkosViewWrite(Mat, MatScalarKokkosView *); 2289d13fa56SJunchao Zhang PETSC_INTERN PetscErrorCode MatInvertVariableBlockDiagonal_SeqAIJKokkos(Mat, const PetscIntKokkosView &, const PetscIntKokkosView &, const PetscIntKokkosView &, PetscScalarKokkosView &, PetscScalarKokkosView &); 229