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