Amesos2 - Direct Sparse Solver Interfaces Version of the Day
Amesos2_MatrixAdapter_decl.hpp
1// @HEADER
2// *****************************************************************************
3// Amesos2: Templated Direct Sparse Solver Package
4//
5// Copyright 2011 NTESS and the Amesos2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10
11#ifndef AMESOS2_MATRIXADAPTER_DECL_HPP
12#define AMESOS2_MATRIXADAPTER_DECL_HPP
13
14#include "Amesos2_config.h"
15
16#include <Teuchos_Comm.hpp>
17#include <Teuchos_ArrayView.hpp>
18#include <Teuchos_VerbosityLevel.hpp>
19#include <Teuchos_FancyOStream.hpp>
20
21#include <Tpetra_ConfigDefs.hpp> // for global_size_t
22
23// #include "Amesos2_ConcreteMatrixAdapter_decl.hpp"
24#include "Amesos2_Util.hpp"
25#include "Amesos2_MatrixTraits.hpp"
26
27namespace Amesos2 {
28
29 template <class M> class ConcreteMatrixAdapter;
30
41 template < class Matrix >
43
44 public:
45
46 typedef typename MatrixTraits<Matrix>::scalar_t scalar_t;
47 typedef typename MatrixTraits<Matrix>::local_ordinal_t local_ordinal_t;
48 typedef typename MatrixTraits<Matrix>::global_ordinal_t global_ordinal_t;
49 typedef typename MatrixTraits<Matrix>::node_t node_t;
50 typedef Tpetra::global_size_t global_size_t;
51
52 typedef Matrix matrix_t;
53 typedef MatrixAdapter<Matrix> type;
54 typedef ConcreteMatrixAdapter<Matrix> adapter_t;
55 typedef Tpetra::Map<local_ordinal_t, global_ordinal_t, node_t> map_t;
56
57 typedef typename MatrixTraits<Matrix>::global_host_idx_type global_host_idx_t;
58 typedef typename MatrixTraits<Matrix>::global_host_val_type global_host_val_t;
59
60 // template<typename S, typename GO, typename GS, typename Op>
61 // friend class Util::get_cxs_helper<MatrixAdapter<Matrix>,S,GO,GS,Op>;
62 // template<class M, typename S, typename GO, typename GS, typename Op>
63 // friend class Util::get_cxs_helper;
64
65 MatrixAdapter(Teuchos::RCP<Matrix> m);
66
67
99 template<typename KV_S, typename KV_GO, typename KV_GS>
100 void getCrs_kokkos_view(KV_S & nzval,
101 KV_GO & colind,
102 KV_GS & rowptr,
103 global_size_t& nnz,
104 const Teuchos::Ptr<const map_t> rowmap,
105 EStorage_Ordering ordering=ARBITRARY,
106 EDistribution distribution=ROOTED) const; // This was placed as last argument to preserve API
107
108
114 template<typename KV_S, typename KV_GO, typename KV_GS>
115 void getCrs_kokkos_view(KV_S & nzval,
116 KV_GO & colind,
117 KV_GS & rowptr,
118 global_size_t& nnz,
119 EDistribution distribution=ROOTED,
120 EStorage_Ordering ordering=ARBITRARY) const; // This was placed as last argument to preserve API
121
150 template<typename KV_S, typename KV_GO, typename KV_GS>
151 void getCcs_kokkos_view(KV_S & nzval,
152 KV_GO & rowind,
153 KV_GS & colptr,
154 global_size_t& nnz,
155 const Teuchos::Ptr<const map_t> colmap,
156 EStorage_Ordering ordering=ARBITRARY,
157 EDistribution distribution=ROOTED) const; // This was placed as last argument to preserve API
158
164 template<typename KV_S, typename KV_GO, typename KV_GS>
165 void getCcs_kokkos_view(KV_S & nzval,
166 KV_GO & rowind,
167 KV_GS & colptr,
168 global_size_t& nnz,
169 EDistribution distribution=ROOTED,
170 EStorage_Ordering ordering=ARBITRARY) const; // This was placed as last argument to preserve API
171
172
174 const Teuchos::RCP<const Teuchos::Comm<int> > getComm() const
175 {
176 return comm_;
177 }
178
180 global_size_t getGlobalNumRows() const;
181
183 global_size_t getGlobalNumCols() const;
184
186 global_size_t getRowIndexBase() const;
187
189 global_size_t getColumnIndexBase() const;
190
192 global_size_t getGlobalNNZ() const;
193
195 size_t getLocalNumRows() const;
196
198 size_t getLocalNumCols() const;
199
201 size_t getLocalNNZ() const;
202
203 Teuchos::RCP<const map_t>
204 getMap() const {
205 return static_cast<const adapter_t*>(this)->getMap_impl();
206 }
207
208 Teuchos::RCP<const map_t>
209 getRowMap() const {
210 return row_map_;
211 }
212
213 Teuchos::RCP<const map_t>
214 getColMap() const {
215 return col_map_;
216 }
217
218 Teuchos::RCP<const type> get(const Teuchos::Ptr<const map_t> map, EDistribution distribution = ROOTED) const;
219
222 Teuchos::RCP<const type> reindex(Teuchos::RCP<const map_t> &contigRowMap, Teuchos::RCP<const map_t> &contigColMap, const EPhase current_phase) const;
223
225 template<typename KV_S, typename KV_GO, typename KV_GS, typename host_ordinal_type_array, typename host_scalar_type_array>
226 local_ordinal_t gather(KV_S& nzvals, KV_GO& indices, KV_GS& pointers,
227 host_ordinal_type_array &perm_g2l,
228 host_ordinal_type_array &recvCountRows, host_ordinal_type_array &recvDisplRows,
229 host_ordinal_type_array &recvCounts, host_ordinal_type_array &recvDispls,
230 host_ordinal_type_array &transpose_map, host_scalar_type_array &nzvals_t,
231 bool column_major, EPhase current_phase) const;
232
234 std::string description() const;
235
237 void describe(Teuchos::FancyOStream &out,
238 const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const;
239
241 template<typename KV>
242 void returnRowPtr_kokkos_view(KV & view) const;
243
245 template<typename KV>
246 void returnColInd_kokkos_view(KV & view) const;
247
249 template<typename KV>
250 void returnValues_kokkos_view(KV & view) const;
251
252
253 private:
254 template<typename KV_S, typename KV_GO, typename KV_GS>
255 void help_getCrs_kokkos_view(KV_S & nzval,
256 KV_GO & colind,
257 KV_GS & rowptr,
258 global_size_t& nnz,
259 const Teuchos::Ptr<const Tpetra::Map<local_ordinal_t,global_ordinal_t,node_t> > rowmap,
260 EDistribution distribution,
261 EStorage_Ordering ordering,
262 no_special_impl nsi) const;
263
264 template<typename KV_S, typename KV_GO, typename KV_GS>
265 void do_getCrs_kokkos_view(KV_S & nzval,
266 KV_GO & colind,
267 KV_GS & rowptr,
268 global_size_t& nnz,
269 const Teuchos::Ptr<const Tpetra::Map<local_ordinal_t,global_ordinal_t,node_t> > rowmap,
270 EDistribution distribution,
271 EStorage_Ordering ordering,
272 row_access ra) const;
273
274 template<typename KV_S, typename KV_GO, typename KV_GS>
275 void help_getCcs_kokkos_view(KV_S & nzval,
276 KV_GO & colind,
277 KV_GS & rowptr,
278 global_size_t& nnz,
279 const Teuchos::Ptr<const Tpetra::Map<local_ordinal_t,global_ordinal_t,node_t> > rowmap,
280 EDistribution distribution,
281 EStorage_Ordering ordering,
282 no_special_impl nsi) const;
283
284 template<typename KV_S, typename KV_GO, typename KV_GS>
285 void do_getCcs_kokkos_view(KV_S & nzval,
286 KV_GO & rowind,
287 KV_GS & colptr,
288 global_size_t& nnz,
289 const Teuchos::Ptr<const Tpetra::Map<local_ordinal_t,global_ordinal_t,node_t> > rowmap,
290 EDistribution distribution,
291 EStorage_Ordering ordering,
292 row_access ra) const;
293
294 protected:
295 // These methods will link to concrete implementations, and may
296 // also be used by them
297
304 template<typename KV_GO, typename KV_S>
305 void getGlobalRowCopy_kokkos_view(global_ordinal_t row,
306 KV_GO & indices,
307 KV_S & vals,
308 size_t& nnz) const;
309
310 size_t getMaxRowNNZ() const;
311
312 size_t getMaxColNNZ() const;
313
314 size_t getGlobalRowNNZ(global_ordinal_t row) const;
315
316 size_t getLocalRowNNZ(local_ordinal_t row) const;
317
318 size_t getGlobalColNNZ(global_ordinal_t col) const;
319
320 size_t getLocalColNNZ(local_ordinal_t col) const;
321
322 bool isLocallyIndexed() const;
323
324 bool isGloballyIndexed() const;
325
326 protected:
327 const Teuchos::RCP<const Matrix> mat_;
328
329 mutable Teuchos::RCP<const Tpetra::Map<local_ordinal_t,global_ordinal_t,node_t> > row_map_;
330
331 mutable Teuchos::RCP<const Tpetra::Map<local_ordinal_t,global_ordinal_t,node_t> > col_map_;
332
333 mutable Teuchos::RCP<const Teuchos::Comm<int> > comm_;
334 }; // end class MatrixAdapter
335
336
337 // Factory creation method
338 template <class Matrix>
339 Teuchos::RCP<MatrixAdapter<Matrix> >
340 createMatrixAdapter(Teuchos::RCP<Matrix> m);
341
342 template <class Matrix>
343 Teuchos::RCP<const MatrixAdapter<Matrix> >
344 createConstMatrixAdapter(Teuchos::RCP<const Matrix> m);
345
346} // end namespace Amesos2
347
348#endif // AMESOS2_MATRIXADAPTER_DECL_HPP
EDistribution
Definition Amesos2_TypeDecl.hpp:89
EStorage_Ordering
Definition Amesos2_TypeDecl.hpp:107
Utility functions for Amesos2.
A Matrix adapter interface for Amesos2.
Definition Amesos2_MatrixAdapter_decl.hpp:42
size_t getLocalNumRows() const
Get the number of rows local to the calling process.
Definition Amesos2_MatrixAdapter_def.hpp:145
void getGlobalRowCopy_kokkos_view(global_ordinal_t row, KV_GO &indices, KV_S &vals, size_t &nnz) const
Definition Amesos2_MatrixAdapter_def.hpp:444
const Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the Teuchos::Comm object associated with this matrix.
Definition Amesos2_MatrixAdapter_decl.hpp:174
Teuchos::RCP< const type > reindex(Teuchos::RCP< const map_t > &contigRowMap, Teuchos::RCP< const map_t > &contigColMap, const EPhase current_phase) const
Definition Amesos2_MatrixAdapter_def.hpp:519
local_ordinal_t gather(KV_S &nzvals, KV_GO &indices, KV_GS &pointers, host_ordinal_type_array &perm_g2l, host_ordinal_type_array &recvCountRows, host_ordinal_type_array &recvDisplRows, host_ordinal_type_array &recvCounts, host_ordinal_type_array &recvDispls, host_ordinal_type_array &transpose_map, host_scalar_type_array &nzvals_t, bool column_major, EPhase current_phase) const
Gather matrix to MPI-0.
Definition Amesos2_MatrixAdapter_def.hpp:527
global_size_t getColumnIndexBase() const
Get the indexbase for the column map.
Definition Amesos2_MatrixAdapter_def.hpp:129
std::string description() const
Returns a short description of this Solver.
Definition Amesos2_MatrixAdapter_def.hpp:167
void getCcs_kokkos_view(KV_S &nzval, KV_GO &rowind, KV_GS &colptr, global_size_t &nnz, EDistribution distribution=ROOTED, EStorage_Ordering ordering=ARBITRARY) const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Describes of this matrix adapter with some level of verbosity.
Definition Amesos2_MatrixAdapter_def.hpp:177
size_t getLocalNumCols() const
Get the number of columns local to the calling process.
Definition Amesos2_MatrixAdapter_def.hpp:152
global_size_t getGlobalNumRows() const
Get the number of rows in this matrix.
Definition Amesos2_MatrixAdapter_def.hpp:106
global_size_t getGlobalNNZ() const
Get the global number of non-zeros in this sparse matrix.
Definition Amesos2_MatrixAdapter_def.hpp:138
void returnColInd_kokkos_view(KV &view) const
Return kokkos view of CRS column indices of matrixA_.
Definition Amesos2_MatrixAdapter_def.hpp:193
void getCrs_kokkos_view(KV_S &nzval, KV_GO &colind, KV_GS &rowptr, global_size_t &nnz, const Teuchos::Ptr< const map_t > rowmap, EStorage_Ordering ordering=ARBITRARY, EDistribution distribution=ROOTED) const
Gets a compressed-row storage summary of this.
void getCrs_kokkos_view(KV_S &nzval, KV_GO &colind, KV_GS &rowptr, global_size_t &nnz, EDistribution distribution=ROOTED, EStorage_Ordering ordering=ARBITRARY) const
size_t getLocalNNZ() const
Get the local number of non-zeros on this processor.
Definition Amesos2_MatrixAdapter_def.hpp:159
global_size_t getGlobalNumCols() const
Get the number of columns in this matrix.
Definition Amesos2_MatrixAdapter_def.hpp:113
void returnValues_kokkos_view(KV &view) const
Return kokkos view of CRS values of matrixA_.
Definition Amesos2_MatrixAdapter_def.hpp:200
void getCcs_kokkos_view(KV_S &nzval, KV_GO &rowind, KV_GS &colptr, global_size_t &nnz, const Teuchos::Ptr< const map_t > colmap, EStorage_Ordering ordering=ARBITRARY, EDistribution distribution=ROOTED) const
Gets a compressed-column storage summary of this.
void returnRowPtr_kokkos_view(KV &view) const
Return kokkos view of CRS row pointer of matrixA_.
Definition Amesos2_MatrixAdapter_def.hpp:186
global_size_t getRowIndexBase() const
Get the indexbase for the row map.
Definition Amesos2_MatrixAdapter_def.hpp:120
EPhase
Used to indicate a phase in the direct solution.
Definition Amesos2_TypeDecl.hpp:31