Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_OverlappingRowMatrix_decl.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_OVERLAPPINGROWMATRIX_DECL_HPP
11#define IFPACK2_OVERLAPPINGROWMATRIX_DECL_HPP
12
13#include "Ifpack2_Details_RowMatrix.hpp"
14#include "Tpetra_CrsMatrix_decl.hpp" // only need the declaration here
15#include "Tpetra_Import_decl.hpp"
16#include "Tpetra_Map_decl.hpp"
17#include <type_traits>
18
19namespace Ifpack2 {
20
24template <class MatrixType>
25class OverlappingRowMatrix : virtual public Ifpack2::Details::RowMatrix<MatrixType> {
26 public:
28
29 typedef typename MatrixType::scalar_type scalar_type;
30 typedef typename MatrixType::local_ordinal_type local_ordinal_type;
31 typedef typename MatrixType::global_ordinal_type global_ordinal_type;
32 typedef typename MatrixType::node_type node_type;
33 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
34 typedef typename MatrixType::global_inds_host_view_type global_inds_host_view_type;
35 typedef typename MatrixType::local_inds_host_view_type local_inds_host_view_type;
36 typedef typename MatrixType::values_host_view_type values_host_view_type;
37
38 typedef typename MatrixType::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
39 typedef typename MatrixType::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
40 typedef typename MatrixType::nonconst_values_host_view_type nonconst_values_host_view_type;
41
42 using row_matrix_type = Tpetra::RowMatrix<scalar_type, local_ordinal_type,
43 global_ordinal_type, node_type>;
44 using crs_matrix_type = Tpetra::CrsMatrix<scalar_type, local_ordinal_type,
45 global_ordinal_type, node_type>;
46
47 // device typedefs
48 typedef typename MatrixType::node_type::device_type device_type;
49 typedef typename device_type::execution_space execution_space;
50 typedef typename MatrixType::local_inds_device_view_type local_inds_device_view_type;
51 typedef typename MatrixType::global_inds_device_view_type global_inds_device_view_type;
52 typedef typename MatrixType::values_device_view_type values_device_view_type;
53
54 static_assert(std::is_same<MatrixType, row_matrix_type>::value, "Ifpack2::OverlappingRowMatrix: The template parameter MatrixType must be a Tpetra::RowMatrix specialization. Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore. The constructor can take either a RowMatrix or a CrsMatrix just fine.");
55
56 typedef typename row_matrix_type::mag_type mag_type;
57
59
61
72 OverlappingRowMatrix(const Teuchos::RCP<const row_matrix_type> &A,
73 const int overlapLevel);
74
77
79
81
83 virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm() const;
84
86 virtual Teuchos::RCP<const Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> >
87 getRowMap() const;
88
90 virtual Teuchos::RCP<const Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> >
91 getColMap() const;
92
96 virtual Teuchos::RCP<const Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> >
97 getDomainMap() const;
98
102 virtual Teuchos::RCP<const Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> >
103 getRangeMap() const;
104
106 virtual Teuchos::RCP<const Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type> >
107 getGraph() const;
108
110 virtual global_size_t getGlobalNumRows() const;
111
113 virtual global_size_t getGlobalNumCols() const;
114
116 virtual size_t getLocalNumRows() const;
117
123 virtual size_t getLocalNumCols() const;
124
126 virtual global_ordinal_type getIndexBase() const;
127
129 virtual global_size_t getGlobalNumEntries() const;
130
132 virtual size_t getLocalNumEntries() const;
133
143 virtual size_t getNumEntriesInGlobalRow(global_ordinal_type globalRow) const;
144
154 virtual size_t getNumEntriesInLocalRow(local_ordinal_type localRow) const;
155
157 virtual size_t getGlobalMaxNumRowEntries() const;
158
160 virtual size_t getLocalMaxNumRowEntries() const;
161
163 virtual local_ordinal_type getBlockSize() const;
164
166 virtual bool hasColMap() const;
167
169 virtual bool isLocallyIndexed() const;
170
172 virtual bool isGloballyIndexed() const;
173
175 virtual bool isFillComplete() const;
176
178 virtual bool supportsRowViews() const;
179
181
183
185
195 virtual void
196 getGlobalRowCopy(global_ordinal_type GlobalRow,
197 nonconst_global_inds_host_view_type &Indices,
198 nonconst_values_host_view_type &Values,
199 size_t &NumEntries) const;
200
202
212 virtual void
213 getLocalRowCopy(local_ordinal_type LocalRow,
214 nonconst_local_inds_host_view_type &Indices,
215 nonconst_values_host_view_type &Values,
216 size_t &NumEntries) const;
217
219
228 virtual void
229 getGlobalRowView(global_ordinal_type GlobalRow,
230 global_inds_host_view_type &indices,
231 values_host_view_type &values) const;
232
234
243 virtual void
244 getLocalRowView(local_ordinal_type LocalRow,
245 local_inds_host_view_type &indices,
246 values_host_view_type &values) const;
247
249
251 virtual void getLocalDiagCopy(Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &diag) const;
252
254
256
266 virtual void
267 leftScale(const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &x);
268
278 virtual void
279 rightScale(const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &x);
280
282
285 virtual mag_type
286 getFrobeniusNorm() const;
287
289
296 virtual void
297 apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &X,
298 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &Y,
299 Teuchos::ETransp mode = Teuchos::NO_TRANS,
300 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
301 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
302
304 virtual bool hasTransposeApply() const;
305
306 virtual void
307 importMultiVector(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &X,
308 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &OvX,
309 Tpetra::CombineMode CM = Tpetra::INSERT);
310
311 virtual void
312 exportMultiVector(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &OvX,
313 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &X,
314 Tpetra::CombineMode CM = Tpetra::ADD);
315
316 std::string description() const;
317
318 void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const;
319
320 Teuchos::RCP<const crs_matrix_type> getUnderlyingMatrix() const;
321
322 Teuchos::RCP<const crs_matrix_type> getExtMatrix() const;
323
324 Kokkos::View<size_t *, typename OverlappingRowMatrix<MatrixType>::device_type> getExtHaloStarts() const;
325 typename Kokkos::View<size_t *, typename OverlappingRowMatrix<MatrixType>::device_type>::host_mirror_type getExtHaloStartsHost() const;
326
327 void doExtImport();
328
329 private:
330 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
331 typedef Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type> import_type;
332 typedef Tpetra::Export<local_ordinal_type, global_ordinal_type, node_type> export_type;
333 typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type> row_graph_type;
334 typedef Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> vector_type;
335
337 Teuchos::RCP<const crs_matrix_type> A_;
338
339 Tpetra::global_size_t NumGlobalRows_;
340 Tpetra::global_size_t NumGlobalNonzeros_;
341 size_t MaxNumEntries_;
342 int OverlapLevel_;
343
344 // Wrapper matrix objects
345 Teuchos::RCP<const map_type> RowMap_;
346 Teuchos::RCP<const map_type> ColMap_;
347 Teuchos::RCP<const import_type> Importer_;
348
350 Teuchos::RCP<crs_matrix_type> ExtMatrix_;
351 Teuchos::RCP<const map_type> ExtMap_;
352 Teuchos::RCP<const import_type> ExtImporter_;
353 Kokkos::View<size_t *, device_type> ExtHaloStarts_;
354 typename Kokkos::View<size_t *, device_type>::host_mirror_type ExtHaloStarts_h;
355
357 Teuchos::RCP<const row_graph_type> graph_;
359 mutable nonconst_local_inds_host_view_type Indices_;
361 mutable nonconst_values_host_view_type Values_;
362
363}; // class OverlappingRowMatrix
364
365} // namespace Ifpack2
366
367#endif // IFPACK2_OVERLAPPINGROWMATRIX_DECL_HPP
All Ifpack2 implementations of Tpetra::RowMatrix must inherit from this class.
Definition Ifpack2_Details_RowMatrix.hpp:33
Sparse matrix (Tpetra::RowMatrix subclass) with ghost rows.
Definition Ifpack2_OverlappingRowMatrix_decl.hpp:25
virtual void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Computes the operator-multivector application.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:452
virtual bool hasTransposeApply() const
Whether this operator's apply() method can apply the adjoint (transpose).
Definition Ifpack2_OverlappingRowMatrix_def.hpp:501
virtual Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const
The Map that describes the range of this matrix.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:238
virtual void getLocalRowCopy(local_ordinal_type LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the graph. Put into storage allocated by callin...
Definition Ifpack2_OverlappingRowMatrix_def.hpp:355
virtual bool isLocallyIndexed() const
Whether this matrix is locally indexed.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:330
virtual size_t getNumEntriesInGlobalRow(global_ordinal_type globalRow) const
The number of entries in the given global row that are owned by the calling process.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:287
virtual void getLocalRowView(local_ordinal_type LocalRow, local_inds_host_view_type &indices, values_host_view_type &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:389
virtual void getGlobalRowCopy(global_ordinal_type GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:346
virtual Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getColMap() const
The Map that describes the distribution of columns over processes.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:217
~OverlappingRowMatrix()=default
Destructor.
virtual size_t getGlobalMaxNumRowEntries() const
The maximum number of entries in any row on any process.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:310
virtual global_size_t getGlobalNumCols() const
The global number of columns in this matrix.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:254
virtual bool isGloballyIndexed() const
Whether this matrix is globally indexed.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:335
virtual global_ordinal_type getIndexBase() const
The index base for global indices for this matrix.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:270
virtual size_t getLocalNumEntries() const
The number of entries in this matrix owned by the calling process.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:280
virtual Teuchos::RCP< const Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > > getGraph() const
This matrix's graph.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:244
virtual size_t getLocalNumCols() const
The number of columns owned by the calling process.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:264
virtual size_t getLocalNumRows() const
The number of rows owned by the calling process.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:259
virtual void getGlobalRowView(global_ordinal_type GlobalRow, global_inds_host_view_type &indices, values_host_view_type &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:371
virtual size_t getNumEntriesInLocalRow(local_ordinal_type localRow) const
The number of entries in the given local row that are owned by the calling process.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:299
virtual void rightScale(const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &x)
Scales the RowMatrix on the right with the Vector x.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:440
virtual local_ordinal_type getBlockSize() const
The number of degrees of freedom per mesh point.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:320
virtual bool hasColMap() const
Whether this matrix has a column Map.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:325
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the matrix is distributed.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:204
virtual global_size_t getGlobalNumRows() const
The global number of rows in this matrix.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:249
virtual bool supportsRowViews() const
true if row views are supported, else false.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:506
virtual global_size_t getGlobalNumEntries() const
The global number of entries in this matrix.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:275
virtual Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRowMap() const
The Map that describes the distribution of rows over processes.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:210
virtual void leftScale(const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &x)
Scales the RowMatrix on the left with the Vector x.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:434
virtual size_t getLocalMaxNumRowEntries() const
The maximum number of entries in any row on the calling process.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:315
virtual bool isFillComplete() const
true if fillComplete() has been called, else false.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:340
virtual mag_type getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:446
virtual void getLocalDiagCopy(Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:404
virtual Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const
The Map that describes the domain of this matrix.
Definition Ifpack2_OverlappingRowMatrix_def.hpp:224
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40