Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_LocalFilter_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_LOCALFILTER_DECL_HPP
11#define IFPACK2_LOCALFILTER_DECL_HPP
12
13#include "Ifpack2_ConfigDefs.hpp"
14#include "Ifpack2_Details_RowMatrix.hpp"
15#include "Tpetra_CrsGraph.hpp"
16#include <type_traits>
17#include <vector>
18
19namespace Ifpack2 {
20
126template <class MatrixType>
127class LocalFilter : virtual public Ifpack2::Details::RowMatrix<MatrixType>,
128 virtual public Teuchos::Describable {
129 private:
130 // Tpetra needs C++11 now because Kokkos needs C++11 now.
131 // Thus, Ifpack2 needs C++11.
132 static_assert(std::is_same<
133 MatrixType,
134 Tpetra::RowMatrix<
135 typename MatrixType::scalar_type,
136 typename MatrixType::local_ordinal_type,
137 typename MatrixType::global_ordinal_type,
138 typename MatrixType::node_type> >::value,
139 "Ifpack2::LocalFilter: MatrixType must be a Tpetra::RowMatrix specialization.");
140
141 public:
143
144
146 typedef typename MatrixType::scalar_type scalar_type;
147
149 typedef typename MatrixType::local_ordinal_type local_ordinal_type;
150
152 typedef typename MatrixType::global_ordinal_type global_ordinal_type;
153
155 typedef typename MatrixType::node_type node_type;
156
157 typedef typename MatrixType::global_inds_host_view_type global_inds_host_view_type;
158 typedef typename MatrixType::local_inds_host_view_type local_inds_host_view_type;
159 typedef typename MatrixType::values_host_view_type values_host_view_type;
160
161 typedef typename MatrixType::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
162 typedef typename MatrixType::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
163 typedef typename MatrixType::nonconst_values_host_view_type nonconst_values_host_view_type;
164
166 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
167
169 typedef Tpetra::RowMatrix<scalar_type,
172 node_type>
174
176 typedef Tpetra::RowGraph<local_ordinal_type,
178 node_type>
180
182 typedef Tpetra::Map<local_ordinal_type,
184 node_type>
186
187 typedef typename row_matrix_type::mag_type mag_type;
188
190
192
194 virtual std::string description() const;
195
197 virtual void
198 describe(Teuchos::FancyOStream &out,
199 const Teuchos::EVerbosityLevel verbLevel =
200 Teuchos::Describable::verbLevel_default) const;
201
203
205
211 explicit LocalFilter(const Teuchos::RCP<const row_matrix_type> &A);
212
214 virtual ~LocalFilter();
215
217
219
221 virtual Teuchos::RCP<const Teuchos::Comm<int> > getComm() const;
222
224 virtual Teuchos::RCP<const map_type> getRowMap() const;
225
227 virtual Teuchos::RCP<const map_type> getColMap() const;
228
230 virtual Teuchos::RCP<const map_type> getDomainMap() const;
231
233 virtual Teuchos::RCP<const map_type> getRangeMap() const;
234
236 virtual Teuchos::RCP<const Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type> >
237 getGraph() const;
238
240 virtual global_size_t getGlobalNumRows() const;
241
243 virtual global_size_t getGlobalNumCols() const;
244
246 virtual size_t getLocalNumRows() const;
247
249 virtual size_t getLocalNumCols() const;
250
252 virtual global_ordinal_type getIndexBase() const;
253
255 virtual global_size_t getGlobalNumEntries() const;
256
258 virtual size_t getLocalNumEntries() const;
259
261 virtual local_ordinal_type getBlockSize() const;
262
268 virtual size_t getNumEntriesInGlobalRow(global_ordinal_type globalRow) const;
269
275 virtual size_t getNumEntriesInLocalRow(local_ordinal_type localRow) const;
276
278 virtual size_t getGlobalMaxNumRowEntries() const;
279
281 virtual size_t getLocalMaxNumRowEntries() const;
282
284 virtual bool hasColMap() const;
285
287 virtual bool isLocallyIndexed() const;
288
290 virtual bool isGloballyIndexed() const;
291
293 virtual bool isFillComplete() const;
294
296 virtual bool supportsRowViews() const;
297
299
301
315 virtual void
317 nonconst_global_inds_host_view_type &Indices,
318 nonconst_values_host_view_type &Values,
319 size_t &NumEntries) const;
320
334 virtual void
336 nonconst_local_inds_host_view_type &Indices,
337 nonconst_values_host_view_type &Values,
338 size_t &NumEntries) const;
339
341
351 virtual void
353 global_inds_host_view_type &indices,
354 values_host_view_type &values) const;
355
357
367 virtual void
369 local_inds_host_view_type &indices,
370 values_host_view_type &values) const;
371
378 virtual void
379 getLocalDiagCopy(Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &diag) const;
380
382
384
394 virtual void leftScale(const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &x);
395
405 virtual void rightScale(const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &x);
406
415 virtual mag_type getFrobeniusNorm() const;
416
423 virtual void
424 apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &X,
425 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &Y,
426 Teuchos::ETransp mode = Teuchos::NO_TRANS,
427 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
428 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
429
431 virtual bool hasTransposeApply() const;
432
434 virtual Teuchos::RCP<const row_matrix_type> getUnderlyingMatrix() const;
435
437 private:
439 typedef Tpetra::CrsGraph<local_ordinal_type,
441 node_type>
442 crs_graph_type;
443
445 void
446 applyNonAliased(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &X,
447 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> &Y,
448 Teuchos::ETransp mode,
449 scalar_type alpha,
450 scalar_type beta) const;
451
464 static bool
465 mapPairIsFitted(const map_type &map1, const map_type &map2);
466
470 // If both pairs of Maps of the original matrix A are fitted on this
471 // process, then this process can use a fast "view" implementation.
472 static bool
473 mapPairsAreFitted(const row_matrix_type &A);
474
476 Teuchos::RCP<const row_matrix_type> A_;
477
479 Teuchos::RCP<const map_type> localRowMap_;
480
482 Teuchos::RCP<const map_type> localDomainMap_;
483
485 Teuchos::RCP<const map_type> localRangeMap_;
486
488 mutable Teuchos::RCP<const row_graph_type> local_graph_;
489
491 size_t NumNonzeros_;
493 size_t MaxNumEntries_;
495 size_t MaxNumEntriesA_;
497 std::vector<size_t> NumEntries_;
498
500 mutable nonconst_local_inds_host_view_type localIndices_;
501 mutable nonconst_local_inds_host_view_type localIndicesForGlobalCopy_;
503 mutable nonconst_values_host_view_type Values_;
504
505}; // class LocalFilter
506
507} // namespace Ifpack2
508
509#endif /* IFPACK2_LOCALFILTER_DECL_HPP */
All Ifpack2 implementations of Tpetra::RowMatrix must inherit from this class.
Definition Ifpack2_Details_RowMatrix.hpp:33
Access only local rows and columns of a sparse matrix.
Definition Ifpack2_LocalFilter_decl.hpp:128
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object to the given output stream.
Definition Ifpack2_LocalFilter_def.hpp:797
virtual Teuchos::RCP< const Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > > getGraph() const
The (locally filtered) matrix's graph.
Definition Ifpack2_LocalFilter_def.hpp:233
virtual Teuchos::RCP< const map_type > getDomainMap() const
Returns the Map that describes the domain distribution in this matrix.
Definition Ifpack2_LocalFilter_def.hpp:217
virtual size_t getGlobalMaxNumRowEntries() const
The maximum number of entries across all rows/columns on all processes.
Definition Ifpack2_LocalFilter_def.hpp:338
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_LocalFilter_def.hpp:547
virtual void getGlobalRowCopy(global_ordinal_type GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const
Get the entries in the given row, using global indices.
Definition Ifpack2_LocalFilter_def.hpp:369
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Type of the Tpetra::RowMatrix specialization that this class uses.
Definition Ifpack2_LocalFilter_decl.hpp:173
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_LocalFilter_decl.hpp:166
virtual bool isGloballyIndexed() const
Whether the underlying sparse matrix is globally (opposite of locally) indexed.
Definition Ifpack2_LocalFilter_def.hpp:358
virtual Teuchos::RCP< const row_matrix_type > getUnderlyingMatrix() const
Return matrix that LocalFilter was built on.
Definition Ifpack2_LocalFilter_def.hpp:836
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_LocalFilter_decl.hpp:146
virtual Teuchos::RCP< const map_type > getRangeMap() const
Returns the Map that describes the range distribution in this matrix.
Definition Ifpack2_LocalFilter_def.hpp:225
virtual mag_type getFrobeniusNorm() const
The Frobenius norm of the (locally filtered) matrix.
Definition Ifpack2_LocalFilter_def.hpp:746
virtual size_t getNumEntriesInGlobalRow(global_ordinal_type globalRow) const
The current number of entries on this node in the specified global row.
Definition Ifpack2_LocalFilter_def.hpp:302
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition Ifpack2_LocalFilter_decl.hpp:155
virtual global_size_t getGlobalNumCols() const
The number of global columns in this matrix.
Definition Ifpack2_LocalFilter_def.hpp:264
virtual bool hasColMap() const
Whether this matrix has a well-defined column Map.
Definition Ifpack2_LocalFilter_def.hpp:348
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_LocalFilter_def.hpp:516
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_LocalFilter_decl.hpp:152
virtual bool isFillComplete() const
Returns true if fillComplete() has been called.
Definition Ifpack2_LocalFilter_def.hpp:363
virtual bool isLocallyIndexed() const
Whether the underlying sparse matrix is locally (opposite of globally) indexed.
Definition Ifpack2_LocalFilter_def.hpp:353
virtual size_t getLocalNumCols() const
The number of columns in the (locally filtered) matrix.
Definition Ifpack2_LocalFilter_def.hpp:274
virtual Teuchos::RCP< const map_type > getColMap() const
Returns the Map that describes the column distribution in this matrix.
Definition Ifpack2_LocalFilter_def.hpp:209
virtual bool hasTransposeApply() const
Whether this operator supports applying the transpose or conjugate transpose.
Definition Ifpack2_LocalFilter_def.hpp:735
virtual void getLocalDiagCopy(Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &diag) const
Get the diagonal entries of the (locally filtered) matrix.
Definition Ifpack2_LocalFilter_def.hpp:534
virtual size_t getLocalMaxNumRowEntries() const
The maximum number of entries across all rows/columns on this process.
Definition Ifpack2_LocalFilter_def.hpp:343
virtual global_ordinal_type getIndexBase() const
Returns the index base for global indices for this matrix.
Definition Ifpack2_LocalFilter_def.hpp:280
virtual bool supportsRowViews() const
Returns true if RowViews are supported.
Definition Ifpack2_LocalFilter_def.hpp:740
virtual std::string description() const
A one-line description of this object.
Definition Ifpack2_LocalFilter_def.hpp:780
virtual global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
Definition Ifpack2_LocalFilter_def.hpp:285
virtual local_ordinal_type getBlockSize() const
The number of degrees of freedom per mesh point.
Definition Ifpack2_LocalFilter_def.hpp:295
virtual ~LocalFilter()
Destructor.
Definition Ifpack2_LocalFilter_def.hpp:189
virtual size_t getLocalNumRows() const
The number of rows owned on the calling process.
Definition Ifpack2_LocalFilter_def.hpp:269
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
Compute Y = beta*Y + alpha*A_local*X.
Definition Ifpack2_LocalFilter_def.hpp:561
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_LocalFilter_decl.hpp:149
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_LocalFilter_def.hpp:554
virtual void getLocalRowCopy(local_ordinal_type LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const
Get the entries in the given row, using local indices.
Definition Ifpack2_LocalFilter_def.hpp:410
virtual size_t getNumEntriesInLocalRow(local_ordinal_type localRow) const
The current number of entries on this node in the specified local row.
Definition Ifpack2_LocalFilter_def.hpp:319
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_LocalFilter_def.hpp:525
virtual global_size_t getGlobalNumRows() const
The number of global rows in this matrix.
Definition Ifpack2_LocalFilter_def.hpp:259
Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > map_type
Type of the Tpetra::Map specialization that this class uses.
Definition Ifpack2_LocalFilter_decl.hpp:185
virtual size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
Definition Ifpack2_LocalFilter_def.hpp:290
Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > row_graph_type
Type of the Tpetra::RowGraph specialization that this class uses.
Definition Ifpack2_LocalFilter_decl.hpp:179
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the communicator.
Definition Ifpack2_LocalFilter_def.hpp:193
virtual Teuchos::RCP< const map_type > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
Definition Ifpack2_LocalFilter_def.hpp:201
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40