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