43 #ifndef IFPACK2_SINGLETONFILTER_DECL_HPP
44 #define IFPACK2_SINGLETONFILTER_DECL_HPP
46 #include "Ifpack2_ConfigDefs.hpp"
47 #include "Tpetra_ConfigDefs.hpp"
48 #include "Ifpack2_Details_RowMatrix.hpp"
49 #include "Tpetra_MultiVector.hpp"
50 #include "Teuchos_RefCountPtr.hpp"
52 #include <type_traits>
63 template<
class MatrixType>
67 typedef typename MatrixType::scalar_type Scalar;
68 typedef typename MatrixType::local_ordinal_type LocalOrdinal;
69 typedef typename MatrixType::global_ordinal_type GlobalOrdinal;
70 typedef typename MatrixType::node_type Node;
71 typedef typename MatrixType::global_inds_host_view_type global_inds_host_view_type;
72 typedef typename MatrixType::local_inds_host_view_type local_inds_host_view_type;
73 typedef typename MatrixType::values_host_view_type values_host_view_type;
75 typedef typename MatrixType::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
76 typedef typename MatrixType::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
77 typedef typename MatrixType::nonconst_values_host_view_type nonconst_values_host_view_type;
80 typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> row_matrix_type;
81 typedef typename row_matrix_type::mag_type mag_type;
83 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
"Ifpack2::SingletonFilter: The template parameter MatrixType must be a Tpetra::RowMatrix specialization. Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore.");
189 nonconst_global_inds_host_view_type &Indices,
190 nonconst_values_host_view_type &Values,
191 size_t& NumEntries)
const;
206 nonconst_local_inds_host_view_type &Indices,
207 nonconst_values_host_view_type &Values,
208 size_t& NumEntries)
const;
222 global_inds_host_view_type &indices,
223 values_host_view_type &values)
const;
237 local_inds_host_view_type & indices,
238 values_host_view_type & values)
const;
243 virtual void getLocalDiagCopy(Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &diag)
const;
259 virtual void leftScale(
const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x);
270 virtual void rightScale(
const Tpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& x);
284 virtual void apply(
const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &X,
285 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y,
294 virtual void SolveSingletons(
const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
295 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS);
297 template <
class DomainScalar,
class RangeScalar>
298 void SolveSingletonsTempl(
const Tpetra::MultiVector<DomainScalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
299 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS);
302 virtual void CreateReducedRHS(
const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS,
303 const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
304 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS);
306 template <
class DomainScalar,
class RangeScalar>
307 void CreateReducedRHSTempl(
const Tpetra::MultiVector<DomainScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS,
308 const Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
309 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS);
312 virtual void UpdateLHS(
const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedLHS,
313 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS);
315 template <
class DomainScalar,
class RangeScalar>
316 void UpdateLHSTempl(
const Tpetra::MultiVector<DomainScalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedLHS,
317 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS);
330 size_t NumSingletons_;
332 std::vector<LocalOrdinal> SingletonIndex_;
334 std::vector<LocalOrdinal> Reorder_;
336 std::vector<LocalOrdinal> InvReorder_;
342 size_t MaxNumEntries_;
344 size_t MaxNumEntriesA_;
346 std::vector<size_t> NumEntries_;
348 mutable nonconst_local_inds_host_view_type Indices_;
350 mutable nonconst_values_host_view_type Values_;
virtual bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Definition: Ifpack2_SingletonFilter_def.hpp:414
virtual mag_type getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:514
virtual void rightScale(const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Scales the RowMatrix on the right with the Vector x.
Definition: Ifpack2_SingletonFilter_def.hpp:363
virtual ~SingletonFilter()
Destructor.
Definition: Ifpack2_SingletonFilter_def.hpp:143
virtual Teuchos::RCP< const Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const
Returns the RowGraph associated with this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:193
virtual void getLocalRowCopy(LocalOrdinal 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_SingletonFilter_def.hpp:307
All Ifpack2 implementations of Tpetra::RowMatrix must inherit from this class.
Definition: Ifpack2_Details_RowMatrix.hpp:63
virtual void UpdateLHS(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &ReducedLHS, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &LHS)
Updates a full LHS from a reduces LHS.
Definition: Ifpack2_SingletonFilter_def.hpp:493
virtual LocalOrdinal getBlockSize() const
The number of degrees of freedom per mesh point.
Definition: Ifpack2_SingletonFilter_def.hpp:265
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:157
virtual void CreateReducedRHS(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &LHS, const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &ReducedRHS)
Creates a RHS for the reduced singleton-free system.
Definition: Ifpack2_SingletonFilter_def.hpp:455
virtual size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
Definition: Ifpack2_SingletonFilter_def.hpp:259
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map that describes the domain distribution in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:175
virtual bool supportsRowViews() const
Returns true if RowViews are supported.
Definition: Ifpack2_SingletonFilter_def.hpp:420
virtual global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:199
virtual void getGlobalRowCopy(GlobalOrdinal 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_SingletonFilter_def.hpp:296
virtual GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:223
virtual bool isFillComplete() const
Returns true if fillComplete() has been called.
Definition: Ifpack2_SingletonFilter_def.hpp:289
virtual global_size_t getGlobalNumCols() const
Returns the number of global columns in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:205
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map that describes the range distribution in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:184
virtual bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Definition: Ifpack2_SingletonFilter_def.hpp:283
virtual global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:229
virtual bool hasColMap() const
Indicates whether this matrix has a well-defined column map.
Definition: Ifpack2_SingletonFilter_def.hpp:271
virtual bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
Definition: Ifpack2_SingletonFilter_def.hpp:277
virtual void getLocalRowView(LocalOrdinal 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_SingletonFilter_def.hpp:341
virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
Definition: Ifpack2_SingletonFilter_def.hpp:247
virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries on this node in the specified global row.
Definition: Ifpack2_SingletonFilter_def.hpp:241
virtual size_t getLocalNumCols() const
Returns the number of columns needed to apply the forward operator on this node, i.e., the number of elements listed in the column map.
Definition: Ifpack2_SingletonFilter_def.hpp:217
virtual size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:235
virtual size_t getLocalNumRows() const
Returns the number of rows owned on the calling node.
Definition: Ifpack2_SingletonFilter_def.hpp:211
virtual size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
Definition: Ifpack2_SingletonFilter_def.hpp:253
virtual void getLocalDiagCopy(Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
Definition: Ifpack2_SingletonFilter_def.hpp:350
virtual void getGlobalRowView(GlobalOrdinal 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_SingletonFilter_def.hpp:332
SingletonFilter(const Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Matrix)
Constructor.
Definition: Ifpack2_SingletonFilter_def.hpp:56
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the communicator.
Definition: Ifpack2_SingletonFilter_def.hpp:147
virtual void apply(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
Computes the operator-multivector application.
Definition: Ifpack2_SingletonFilter_def.hpp:369
virtual void SolveSingletons(const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &RHS, Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &LHS)
Solve the singleton components of the linear system.
Definition: Ifpack2_SingletonFilter_def.hpp:426
virtual Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:166
Filter based on matrix entries.
Definition: Ifpack2_SingletonFilter_decl.hpp:64
virtual void leftScale(const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Scales the RowMatrix on the left with the Vector x.
Definition: Ifpack2_SingletonFilter_def.hpp:357