xref: /petsc/src/mat/impls/aij/seq/kokkos/aijkok.hpp (revision d2522c19e8fa9bca20aaca277941d9a63e71db6a)
1 #ifndef __SEQAIJKOKKOSIMPL_HPP
2 #define __SEQAIJKOKKOSIMPL_HPP
3 
4 #include <petscaijdevice.h>
5 #include <petsc/private/vecimpl_kokkos.hpp>
6 #include <../src/mat/impls/aij/seq/aij.h>
7 #include <KokkosSparse_CrsMatrix.hpp>
8 #include <KokkosSparse_spiluk.hpp>
9 
10 /*
11    Kokkos::View<struct _n_SplitCSRMat,DefaultMemorySpace> is not handled correctly so we define SplitCSRMat
12    for the singular purpose of working around this.
13 */
14 typedef struct _n_SplitCSRMat SplitCSRMat;
15 
16 using MatRowMapType = PetscInt;
17 using MatColIdxType = PetscInt;
18 using MatScalarType = PetscScalar;
19 
20 template <class MemorySpace>
21 using KokkosCsrMatrixType = typename KokkosSparse::CrsMatrix<MatScalarType, MatColIdxType, MemorySpace, void /* MemoryTraits */, MatRowMapType>;
22 template <class MemorySpace>
23 using KokkosCsrGraphType = typename KokkosCsrMatrixType<MemorySpace>::staticcrsgraph_type;
24 
25 using KokkosCsrGraph     = KokkosCsrGraphType<DefaultMemorySpace>;
26 using KokkosCsrGraphHost = KokkosCsrGraphType<Kokkos::HostSpace>;
27 
28 using KokkosCsrMatrix     = KokkosCsrMatrixType<DefaultMemorySpace>;
29 using KokkosCsrMatrixHost = KokkosCsrMatrixType<Kokkos::HostSpace>;
30 
31 using MatRowMapKokkosView = KokkosCsrGraph::row_map_type::non_const_type;
32 using MatColIdxKokkosView = KokkosCsrGraph::entries_type::non_const_type;
33 using MatScalarKokkosView = KokkosCsrMatrix::values_type::non_const_type;
34 
35 using MatRowMapKokkosViewHost = KokkosCsrGraphHost::row_map_type::non_const_type;
36 using MatColIdxKokkosViewHost = KokkosCsrGraphHost::entries_type::non_const_type;
37 using MatScalarKokkosViewHost = KokkosCsrMatrixHost::values_type::non_const_type;
38 
39 using ConstMatRowMapKokkosView = KokkosCsrGraph::row_map_type::const_type;
40 using ConstMatColIdxKokkosView = KokkosCsrGraph::entries_type::const_type;
41 using ConstMatScalarKokkosView = KokkosCsrMatrix::values_type::const_type;
42 
43 using ConstMatRowMapKokkosViewHost = KokkosCsrGraphHost::row_map_type::const_type;
44 using ConstMatColIdxKokkosViewHost = KokkosCsrGraphHost::entries_type::const_type;
45 using ConstMatScalarKokkosViewHost = KokkosCsrMatrixHost::values_type::const_type;
46 
47 using MatRowMapKokkosDualView = Kokkos::DualView<MatRowMapType *>;
48 using MatColIdxKokkosDualView = Kokkos::DualView<MatColIdxType *>;
49 using MatScalarKokkosDualView = Kokkos::DualView<MatScalarType *>;
50 
51 using KernelHandle = KokkosKernels::Experimental::KokkosKernelsHandle<MatRowMapType, MatColIdxType, MatScalarType, DefaultExecutionSpace, DefaultMemorySpace, DefaultMemorySpace>;
52 
53 using KokkosTeamMemberType = Kokkos::TeamPolicy<DefaultExecutionSpace>::member_type;
54 
55 /* For mat->spptr of a factorized matrix */
56 struct Mat_SeqAIJKokkosTriFactors {
57   MatRowMapKokkosView   iL_d, iU_d, iLt_d, iUt_d;  /* rowmap for L, U, L^t, U^t of A=LU */
58   MatColIdxKokkosView   jL_d, jU_d, jLt_d, jUt_d;  /* column ids */
59   MatScalarKokkosView   aL_d, aU_d, aLt_d, aUt_d;  /* matrix values */
60   KernelHandle          kh, khL, khU, khLt, khUt;  /* Kernel handles for A, L, U, L^t, U^t */
61   PetscBool             transpose_updated;         /* Are L^T, U^T updated wrt L, U*/
62   PetscBool             sptrsv_symbolic_completed; /* Have we completed the symbolic solve for L and U */
63   PetscScalarKokkosView workVector;
64 
65   Mat_SeqAIJKokkosTriFactors(PetscInt n) : transpose_updated(PETSC_FALSE), sptrsv_symbolic_completed(PETSC_FALSE), workVector("workVector", n) { }
66 
67   ~Mat_SeqAIJKokkosTriFactors() { Destroy(); }
68 
69   void Destroy() {
70     kh.destroy_spiluk_handle();
71     khL.destroy_sptrsv_handle();
72     khU.destroy_sptrsv_handle();
73     khLt.destroy_sptrsv_handle();
74     khUt.destroy_sptrsv_handle();
75     transpose_updated = sptrsv_symbolic_completed = PETSC_FALSE;
76   }
77 };
78 
79 /* For mat->spptr of a regular matrix */
80 struct Mat_SeqAIJKokkos {
81   MatRowMapKokkosDualView i_dual;
82   MatColIdxKokkosDualView j_dual;
83   MatScalarKokkosDualView a_dual;
84 
85   MatRowMapKokkosDualView diag_dual; /* Diagonal pointer, built on demand */
86 
87   KokkosCsrMatrix  csrmat;       /* The CSR matrix, used to call KK functions */
88   PetscObjectState nonzerostate; /* State of the nonzero pattern (graph) on device */
89 
90   KokkosCsrMatrix csrmatT, csrmatH;                     /* Transpose and Hermitian of the matrix (built on demand) */
91   PetscBool       transpose_updated, hermitian_updated; /* Are At, Ah updated wrt the matrix? */
92 
93   /* COO stuff */
94   PetscCountKokkosView jmap_d; /* perm[disp+jmap[i]..disp+jmap[i+1]) gives indices of entries in v[] associated with i-th nonzero of the matrix */
95   PetscCountKokkosView perm_d; /* The permutation array in sorting (i,j) by row and then by col */
96 
97   Kokkos::View<PetscInt *>                      i_uncompressed_d;
98   Kokkos::View<PetscInt *>                      colmap_d; // ugh, this is a parallel construct
99   Kokkos::View<SplitCSRMat, DefaultMemorySpace> device_mat_d;
100   Kokkos::View<PetscInt *>                      diag_d; // factorizations
101 
102   /* Construct a nrows by ncols matrix with nnz nonzeros from the given (i,j,a) on host. Caller also specifies a nonzero state */
103   Mat_SeqAIJKokkos(PetscInt nrows, PetscInt ncols, PetscInt nnz, const MatRowMapType *i, MatColIdxType *j, MatScalarType *a, PetscObjectState nzstate, PetscBool copyValues = PETSC_TRUE) {
104     MatScalarKokkosViewHost a_h(a, nnz);
105     MatRowMapKokkosViewHost i_h(const_cast<MatRowMapType *>(i), nrows + 1);
106     MatColIdxKokkosViewHost j_h(j, nnz);
107 
108     auto a_d = Kokkos::create_mirror_view(DefaultMemorySpace(), a_h);
109     auto i_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), i_h);
110     auto j_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), j_h);
111 
112     a_dual = MatScalarKokkosDualView(a_d, a_h);
113     i_dual = MatRowMapKokkosDualView(i_d, i_h);
114     j_dual = MatColIdxKokkosDualView(j_d, j_h);
115 
116     a_dual.modify_host(); /* Since caller provided values on host */
117     if (copyValues) a_dual.sync_device();
118 
119     csrmat            = KokkosCsrMatrix("csrmat", ncols, a_d, KokkosCsrGraph(j_d, i_d));
120     nonzerostate      = nzstate;
121     transpose_updated = hermitian_updated = PETSC_FALSE;
122   }
123 
124   /* Construct with a KokkosCsrMatrix. For performance, only i, j are copied to host, but not the matrix values. */
125   Mat_SeqAIJKokkos(const KokkosCsrMatrix &csr) :
126     csrmat(csr) /* Shallow-copy csr's views to csrmat */
127   {
128     auto                a_d = csr.values;
129     /* Get a non-const version since I don't want to deal with DualView<const T*>, which is not well defined */
130     MatRowMapKokkosView i_d(const_cast<MatRowMapType *>(csr.graph.row_map.data()), csr.graph.row_map.extent(0));
131     auto                j_d = csr.graph.entries;
132     auto                a_h = Kokkos::create_mirror_view(Kokkos::HostSpace(), a_d);
133     auto                i_h = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), i_d);
134     auto                j_h = Kokkos::create_mirror_view_and_copy(Kokkos::HostSpace(), j_d);
135 
136     a_dual = MatScalarKokkosDualView(a_d, a_h);
137     a_dual.modify_device(); /* since we did not copy a_d to a_h, we mark device has the latest data */
138     i_dual = MatRowMapKokkosDualView(i_d, i_h);
139     j_dual = MatColIdxKokkosDualView(j_d, j_h);
140     Init();
141   }
142 
143   Mat_SeqAIJKokkos(PetscInt nrows, PetscInt ncols, PetscInt nnz, MatRowMapKokkosDualView &i, MatColIdxKokkosDualView &j, MatScalarKokkosDualView a) : i_dual(i), j_dual(j), a_dual(a) {
144     csrmat = KokkosCsrMatrix("csrmat", nrows, ncols, nnz, a.view_device(), i.view_device(), j.view_device());
145     Init();
146   }
147 
148   MatScalarType *a_host_data() { return a_dual.view_host().data(); }
149   MatRowMapType *i_host_data() { return i_dual.view_host().data(); }
150   MatColIdxType *j_host_data() { return j_dual.view_host().data(); }
151 
152   MatScalarType *a_device_data() { return a_dual.view_device().data(); }
153   MatRowMapType *i_device_data() { return i_dual.view_device().data(); }
154   MatColIdxType *j_device_data() { return j_dual.view_device().data(); }
155 
156   MatColIdxType nrows() { return csrmat.numRows(); }
157   MatColIdxType ncols() { return csrmat.numCols(); }
158   MatRowMapType nnz() { return csrmat.nnz(); }
159 
160   /* Change the csrmat size to n */
161   void SetColSize(MatColIdxType n) { csrmat = KokkosCsrMatrix("csrmat", n, a_dual.view_device(), csrmat.graph); }
162 
163   void SetUpCOO(const Mat_SeqAIJ *aij) {
164     jmap_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(aij->jmap, aij->nz + 1));
165     perm_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(aij->perm, aij->Atot));
166   }
167 
168   void SetDiagonal(const MatRowMapType *diag) {
169     MatRowMapKokkosViewHost diag_h(const_cast<MatRowMapType *>(diag), nrows());
170     auto                    diag_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), diag_h);
171     diag_dual                      = MatRowMapKokkosDualView(diag_d, diag_h);
172   }
173 
174   /* Shared init stuff */
175   void Init(void) {
176     transpose_updated = hermitian_updated = PETSC_FALSE;
177     nonzerostate                          = 0;
178   }
179 
180   PetscErrorCode DestroyMatTranspose(void) {
181     PetscFunctionBegin;
182     csrmatT = KokkosCsrMatrix(); /* Overwrite with empty matrices */
183     csrmatH = KokkosCsrMatrix();
184     PetscFunctionReturn(0);
185   }
186 };
187 
188 struct MatProductData_SeqAIJKokkos {
189   KernelHandle kh;
190   PetscBool    reusesym;
191   MatProductData_SeqAIJKokkos() : reusesym(PETSC_FALSE) { }
192 };
193 
194 PETSC_INTERN PetscErrorCode MatSetSeqAIJKokkosWithCSRMatrix(Mat, Mat_SeqAIJKokkos *);
195 PETSC_INTERN PetscErrorCode MatCreateSeqAIJKokkosWithCSRMatrix(MPI_Comm, Mat_SeqAIJKokkos *, Mat *);
196 PETSC_INTERN PetscErrorCode MatSeqAIJKokkosMergeMats(Mat, Mat, MatReuse, Mat *);
197 PETSC_INTERN PetscErrorCode MatSeqAIJKokkosSyncDevice(Mat);
198 PETSC_INTERN PetscErrorCode PrintCsrMatrix(const KokkosCsrMatrix &csrmat);
199 PETSC_INTERN PetscErrorCode MatConvert_SeqAIJ_SeqAIJKokkos(Mat, MatType, MatReuse, Mat *);
200 PETSC_INTERN PetscErrorCode MatSeqAIJKokkosModifyDevice(Mat);
201 
202 PETSC_INTERN PetscErrorCode MatSeqAIJGetKokkosView(Mat, MatScalarKokkosView *);
203 PETSC_INTERN PetscErrorCode MatSeqAIJRestoreKokkosView(Mat, MatScalarKokkosView *);
204 PETSC_INTERN PetscErrorCode MatSeqAIJGetKokkosViewWrite(Mat, MatScalarKokkosView *);
205 PETSC_INTERN PetscErrorCode MatSeqAIJRestoreKokkosViewWrite(Mat, MatScalarKokkosView *);
206 #endif
207