xref: /petsc/src/mat/impls/aij/seq/kokkos/aijkok.hpp (revision 9d13fa56c5c6523e02c36edc0e4e22bf2d0334a8)
142550becSJunchao Zhang #ifndef __SEQAIJKOKKOSIMPL_HPP
242550becSJunchao Zhang #define __SEQAIJKOKKOSIMPL_HPP
342550becSJunchao Zhang 
442550becSJunchao Zhang #include <petsc/private/vecimpl_kokkos.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>;
2842550becSJunchao Zhang using KokkosCsrGraphHost = KokkosCsrGraphType<Kokkos::HostSpace>;
2942550becSJunchao Zhang 
3042550becSJunchao Zhang using KokkosCsrMatrix     = KokkosCsrMatrixType<DefaultMemorySpace>;
3142550becSJunchao Zhang using KokkosCsrMatrixHost = KokkosCsrMatrixType<Kokkos::HostSpace>;
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 */
6242550becSJunchao Zhang   KernelHandle          kh, khL, khU, khLt, khUt;  /* Kernel handles for A, L, U, L^t, U^t */
6342550becSJunchao Zhang   PetscBool             transpose_updated;         /* Are L^T, U^T updated wrt L, U*/
6442550becSJunchao Zhang   PetscBool             sptrsv_symbolic_completed; /* Have we completed the symbolic solve for L and U */
6542550becSJunchao Zhang   PetscScalarKokkosView workVector;
6642550becSJunchao Zhang 
679371c9d4SSatish Balay   Mat_SeqAIJKokkosTriFactors(PetscInt n) : transpose_updated(PETSC_FALSE), sptrsv_symbolic_completed(PETSC_FALSE), workVector("workVector", n) { }
6842550becSJunchao Zhang 
6942550becSJunchao Zhang   ~Mat_SeqAIJKokkosTriFactors() { Destroy(); }
7042550becSJunchao Zhang 
71d71ae5a4SJacob Faibussowitsch   void Destroy()
72d71ae5a4SJacob Faibussowitsch   {
7342550becSJunchao Zhang     kh.destroy_spiluk_handle();
7442550becSJunchao Zhang     khL.destroy_sptrsv_handle();
7542550becSJunchao Zhang     khU.destroy_sptrsv_handle();
7642550becSJunchao Zhang     khLt.destroy_sptrsv_handle();
7742550becSJunchao Zhang     khUt.destroy_sptrsv_handle();
7842550becSJunchao Zhang     transpose_updated = sptrsv_symbolic_completed = PETSC_FALSE;
7942550becSJunchao Zhang   }
8042550becSJunchao Zhang };
8142550becSJunchao Zhang 
8242550becSJunchao Zhang /* For mat->spptr of a regular matrix */
8342550becSJunchao Zhang struct Mat_SeqAIJKokkos {
8442550becSJunchao Zhang   MatRowMapKokkosDualView i_dual;
8542550becSJunchao Zhang   MatColIdxKokkosDualView j_dual;
8642550becSJunchao Zhang   MatScalarKokkosDualView a_dual;
8742550becSJunchao Zhang 
88f78ce678SMark Adams   MatRowMapKokkosDualView diag_dual; /* Diagonal pointer, built on demand */
89f78ce678SMark Adams 
9042550becSJunchao Zhang   KokkosCsrMatrix  csrmat;       /* The CSR matrix, used to call KK functions */
9142550becSJunchao Zhang   PetscObjectState nonzerostate; /* State of the nonzero pattern (graph) on device */
9242550becSJunchao Zhang 
9342550becSJunchao Zhang   KokkosCsrMatrix     csrmatT, csrmatH;                     /* Transpose and Hermitian of the matrix (built on demand) */
9442550becSJunchao Zhang   PetscBool           transpose_updated, hermitian_updated; /* Are At, Ah updated wrt the matrix? */
950e3ece09SJunchao Zhang   MatRowMapKokkosView transpose_perm;                       // A permutation array making Ta(i) = Aa(perm(i)), where T = A^t
9642550becSJunchao Zhang 
9742550becSJunchao Zhang   /* Construct a nrows by ncols matrix with nnz nonzeros from the given (i,j,a) on host. Caller also specifies a nonzero state */
98d71ae5a4SJacob Faibussowitsch   Mat_SeqAIJKokkos(PetscInt nrows, PetscInt ncols, PetscInt nnz, const MatRowMapType *i, MatColIdxType *j, MatScalarType *a, PetscObjectState nzstate, PetscBool copyValues = PETSC_TRUE)
99d71ae5a4SJacob Faibussowitsch   {
10042550becSJunchao Zhang     MatScalarKokkosViewHost a_h(a, nnz);
10142550becSJunchao Zhang     MatRowMapKokkosViewHost i_h(const_cast<MatRowMapType *>(i), nrows + 1);
10242550becSJunchao Zhang     MatColIdxKokkosViewHost j_h(j, nnz);
10342550becSJunchao Zhang 
10442550becSJunchao Zhang     auto a_d = Kokkos::create_mirror_view(DefaultMemorySpace(), a_h);
10542550becSJunchao Zhang     auto i_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), i_h);
10642550becSJunchao Zhang     auto j_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), j_h);
10742550becSJunchao Zhang 
10842550becSJunchao Zhang     a_dual = MatScalarKokkosDualView(a_d, a_h);
10942550becSJunchao Zhang     i_dual = MatRowMapKokkosDualView(i_d, i_h);
11042550becSJunchao Zhang     j_dual = MatColIdxKokkosDualView(j_d, j_h);
11142550becSJunchao Zhang 
11242550becSJunchao Zhang     a_dual.modify_host(); /* Since caller provided values on host */
11342550becSJunchao Zhang     if (copyValues) a_dual.sync_device();
11442550becSJunchao Zhang 
11542550becSJunchao Zhang     csrmat            = KokkosCsrMatrix("csrmat", ncols, a_d, KokkosCsrGraph(j_d, i_d));
11642550becSJunchao Zhang     nonzerostate      = nzstate;
11742550becSJunchao Zhang     transpose_updated = hermitian_updated = PETSC_FALSE;
11842550becSJunchao Zhang   }
11942550becSJunchao Zhang 
12042550becSJunchao Zhang   /* Construct with a KokkosCsrMatrix. For performance, only i, j are copied to host, but not the matrix values. */
121d71ae5a4SJacob Faibussowitsch   Mat_SeqAIJKokkos(const KokkosCsrMatrix &csr) : csrmat(csr) /* Shallow-copy csr's views to csrmat */
12242550becSJunchao Zhang   {
12342550becSJunchao Zhang     auto a_d = csr.values;
12442550becSJunchao Zhang     /* Get a non-const version since I don't want to deal with DualView<const T*>, which is not well defined */
12542550becSJunchao Zhang     MatRowMapKokkosView i_d(const_cast<MatRowMapType *>(csr.graph.row_map.data()), csr.graph.row_map.extent(0));
12642550becSJunchao Zhang     auto                j_d = csr.graph.entries;
12742550becSJunchao Zhang     auto                a_h = Kokkos::create_mirror_view(Kokkos::HostSpace(), a_d);
12842550becSJunchao Zhang     auto                i_h = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), i_d);
12942550becSJunchao Zhang     auto                j_h = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), j_d);
13042550becSJunchao Zhang 
13142550becSJunchao Zhang     a_dual = MatScalarKokkosDualView(a_d, a_h);
13242550becSJunchao Zhang     a_dual.modify_device(); /* since we did not copy a_d to a_h, we mark device has the latest data */
13342550becSJunchao Zhang     i_dual = MatRowMapKokkosDualView(i_d, i_h);
13442550becSJunchao Zhang     j_dual = MatColIdxKokkosDualView(j_d, j_h);
13542550becSJunchao Zhang     Init();
13642550becSJunchao Zhang   }
13742550becSJunchao Zhang 
138d71ae5a4SJacob Faibussowitsch   Mat_SeqAIJKokkos(PetscInt nrows, PetscInt ncols, PetscInt nnz, MatRowMapKokkosDualView &i, MatColIdxKokkosDualView &j, MatScalarKokkosDualView a) : i_dual(i), j_dual(j), a_dual(a)
139d71ae5a4SJacob Faibussowitsch   {
14042550becSJunchao Zhang     csrmat = KokkosCsrMatrix("csrmat", nrows, ncols, nnz, a.view_device(), i.view_device(), j.view_device());
14142550becSJunchao Zhang     Init();
14242550becSJunchao Zhang   }
14342550becSJunchao Zhang 
14442550becSJunchao Zhang   MatScalarType *a_host_data() { return a_dual.view_host().data(); }
14542550becSJunchao Zhang   MatRowMapType *i_host_data() { return i_dual.view_host().data(); }
14642550becSJunchao Zhang   MatColIdxType *j_host_data() { return j_dual.view_host().data(); }
14742550becSJunchao Zhang 
14842550becSJunchao Zhang   MatScalarType *a_device_data() { return a_dual.view_device().data(); }
14942550becSJunchao Zhang   MatRowMapType *i_device_data() { return i_dual.view_device().data(); }
15042550becSJunchao Zhang   MatColIdxType *j_device_data() { return j_dual.view_device().data(); }
15142550becSJunchao Zhang 
15242550becSJunchao Zhang   MatColIdxType nrows() { return csrmat.numRows(); }
15342550becSJunchao Zhang   MatColIdxType ncols() { return csrmat.numCols(); }
15442550becSJunchao Zhang   MatRowMapType nnz() { return csrmat.nnz(); }
15542550becSJunchao Zhang 
15642550becSJunchao Zhang   /* Change the csrmat size to n */
15742550becSJunchao Zhang   void SetColSize(MatColIdxType n) { csrmat = KokkosCsrMatrix("csrmat", n, a_dual.view_device(), csrmat.graph); }
15842550becSJunchao Zhang 
159d71ae5a4SJacob Faibussowitsch   void SetDiagonal(const MatRowMapType *diag)
160d71ae5a4SJacob Faibussowitsch   {
161f78ce678SMark Adams     MatRowMapKokkosViewHost diag_h(const_cast<MatRowMapType *>(diag), nrows());
162f78ce678SMark Adams     auto                    diag_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), diag_h);
163f78ce678SMark Adams     diag_dual                      = MatRowMapKokkosDualView(diag_d, diag_h);
164f78ce678SMark Adams   }
165f78ce678SMark Adams 
16642550becSJunchao Zhang   /* Shared init stuff */
167d71ae5a4SJacob Faibussowitsch   void Init(void)
168d71ae5a4SJacob Faibussowitsch   {
16942550becSJunchao Zhang     transpose_updated = hermitian_updated = PETSC_FALSE;
17042550becSJunchao Zhang     nonzerostate                          = 0;
17142550becSJunchao Zhang   }
17242550becSJunchao Zhang 
173d71ae5a4SJacob Faibussowitsch   PetscErrorCode DestroyMatTranspose(void)
174d71ae5a4SJacob Faibussowitsch   {
17542550becSJunchao Zhang     PetscFunctionBegin;
17642550becSJunchao Zhang     csrmatT = KokkosCsrMatrix(); /* Overwrite with empty matrices */
17742550becSJunchao Zhang     csrmatH = KokkosCsrMatrix();
1783ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
17942550becSJunchao Zhang   }
18042550becSJunchao Zhang };
18142550becSJunchao Zhang 
18242550becSJunchao Zhang struct MatProductData_SeqAIJKokkos {
18342550becSJunchao Zhang   KernelHandle kh;
18442550becSJunchao Zhang   PetscBool    reusesym;
18542550becSJunchao Zhang   MatProductData_SeqAIJKokkos() : reusesym(PETSC_FALSE) { }
18642550becSJunchao Zhang };
18742550becSJunchao Zhang 
18842550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatSetSeqAIJKokkosWithCSRMatrix(Mat, Mat_SeqAIJKokkos *);
18942550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatCreateSeqAIJKokkosWithCSRMatrix(MPI_Comm, Mat_SeqAIJKokkos *, Mat *);
19042550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosMergeMats(Mat, Mat, MatReuse, Mat *);
19142550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosSyncDevice(Mat);
1920e3ece09SJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosGetKokkosCsrMatrix(Mat, KokkosCsrMatrix *);
1930e3ece09SJunchao Zhang PETSC_INTERN PetscErrorCode MatCreateSeqAIJKokkosWithKokkosCsrMatrix(MPI_Comm, KokkosCsrMatrix, Mat *);
19442550becSJunchao Zhang PETSC_INTERN PetscErrorCode PrintCsrMatrix(const KokkosCsrMatrix &csrmat);
19542550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJKokkos(Mat, MatType, MatReuse, Mat *);
19642550becSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosModifyDevice(Mat);
1970e3ece09SJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJKokkosGenerateTranspose_Private(Mat, KokkosCsrMatrix *);
198394ed5ebSJunchao Zhang 
199394ed5ebSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJGetKokkosView(Mat, MatScalarKokkosView *);
200394ed5ebSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJRestoreKokkosView(Mat, MatScalarKokkosView *);
201394ed5ebSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJGetKokkosViewWrite(Mat, MatScalarKokkosView *);
202394ed5ebSJunchao Zhang PETSC_INTERN PetscErrorCode MatSeqAIJRestoreKokkosViewWrite(Mat, MatScalarKokkosView *);
203*9d13fa56SJunchao Zhang PETSC_INTERN PetscErrorCode MatInvertVariableBlockDiagonal_SeqAIJKokkos(Mat, const PetscIntKokkosView &, const PetscIntKokkosView &, const PetscIntKokkosView &, PetscScalarKokkosView &, PetscScalarKokkosView &);
204*9d13fa56SJunchao Zhang 
20542550becSJunchao Zhang #endif
206