43#ifndef IFPACK2_SINGLETONFILTER_DEF_HPP
44#define IFPACK2_SINGLETONFILTER_DEF_HPP
45#include "Ifpack2_SingletonFilter_decl.hpp"
47#include "Tpetra_ConfigDefs.hpp"
48#include "Tpetra_RowMatrix.hpp"
49#include "Tpetra_Map.hpp"
50#include "Tpetra_MultiVector.hpp"
51#include "Tpetra_Vector.hpp"
55template<
class MatrixType>
66 if (A_->getComm()->getSize() != 1 || A_->getLocalNumRows() != A_->getGlobalNumRows()) {
67 throw std::runtime_error(
"Ifpack2::SingeltonFilter can be used with Comm().getSize() == 1 only. This class is a tool for Ifpack2_AdditiveSchwarz, and it is not meant to be used otherwise.");
71 size_t NumRowsA_ = A_->getLocalNumRows();
75 MaxNumEntriesA_ = A_->getLocalMaxNumRowEntries();
78 Kokkos::resize(Indices_,MaxNumEntriesA_);
79 Kokkos::resize(Values_,MaxNumEntriesA_);
82 Reorder_.resize(NumRowsA_);
83 Reorder_.assign(Reorder_.size(),-1);
87 for (
size_t i = 0 ; i < NumRowsA_ ; ++i) {
89 A_->getLocalRowCopy(i,Indices_,Values_,Nnz);
91 Reorder_[i] = NumRows_++;
99 InvReorder_.resize(NumRows_);
100 for (
size_t i = 0 ; i < NumRowsA_ ; ++i) {
103 InvReorder_[Reorder_[i]] = i;
105 NumEntries_.resize(NumRows_);
106 SingletonIndex_.resize(NumSingletons_);
111 for (
size_t i = 0 ; i < NumRowsA_ ; ++i) {
113 A_->getLocalRowCopy(i,Indices_,Values_,Nnz);
114 LocalOrdinal ii = Reorder_[i];
116 NumEntries_[ii] = Nnz;
118 if (Nnz > MaxNumEntries_)
119 MaxNumEntries_ = Nnz;
122 SingletonIndex_[count] = i;
128 ReducedMap_ = Teuchos::rcp(
new Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node>(NumRows_,0,A_->getComm()) );
131 Diagonal_ = Teuchos::rcp(
new Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>(ReducedMap_) );
133 Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> DiagonalA(A_->getRowMap());
134 A_->getLocalDiagCopy(DiagonalA);
135 const Teuchos::ArrayRCP<const Scalar> & DiagonalAview = DiagonalA.get1dView();
136 for (
size_t i = 0 ; i < NumRows_ ; ++i) {
137 LocalOrdinal ii = InvReorder_[i];
138 Diagonal_->replaceLocalValue((LocalOrdinal)i,DiagonalAview[ii]);
142template<
class MatrixType>
145template<
class MatrixType>
146Teuchos::RCP<const Teuchos::Comm<int> >
149 return A_->getComm();
153template<
class MatrixType>
154Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
155 typename MatrixType::global_ordinal_type,
156 typename MatrixType::node_type> >
162template<
class MatrixType>
163Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
164 typename MatrixType::global_ordinal_type,
165 typename MatrixType::node_type> >
171template<
class MatrixType>
172Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
173 typename MatrixType::global_ordinal_type,
174 typename MatrixType::node_type> >
180template<
class MatrixType>
181Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
182 typename MatrixType::global_ordinal_type,
183 typename MatrixType::node_type> >
189template<
class MatrixType>
190Teuchos::RCP<
const Tpetra::RowGraph<
typename MatrixType::local_ordinal_type,
191 typename MatrixType::global_ordinal_type,
192 typename MatrixType::node_type> >
195 throw std::runtime_error(
"Ifpack2::SingletonFilter: does not support getGraph.");
198template<
class MatrixType>
204template<
class MatrixType>
210template<
class MatrixType>
216template<
class MatrixType>
222template<
class MatrixType>
225 return A_->getIndexBase();
228template<
class MatrixType>
234template<
class MatrixType>
240template<
class MatrixType>
243 throw std::runtime_error(
"Ifpack2::SingletonFilter does not implement getNumEntriesInGlobalRow.");
246template<
class MatrixType>
249 return NumEntries_[localRow];
252template<
class MatrixType>
255 return MaxNumEntries_;
258template<
class MatrixType>
261 return MaxNumEntries_;
264template<
class MatrixType>
270template<
class MatrixType>
273 return A_->isLocallyIndexed();
276template<
class MatrixType>
279 return A_->isGloballyIndexed();
282template<
class MatrixType>
285 return A_->isFillComplete();
288template<
class MatrixType>
291 nonconst_global_inds_host_view_type &,
292 nonconst_values_host_view_type &,
295 throw std::runtime_error(
"Ifpack2::SingletonFilter does not implement getGlobalRowCopy.");
299template<
class MatrixType>
302 nonconst_local_inds_host_view_type &Indices,
303 nonconst_values_host_view_type &Values,
304 size_t& NumEntries)
const
306 TEUCHOS_TEST_FOR_EXCEPTION((LocalRow < 0 || (
size_t) LocalRow >= NumRows_ || (
size_t) Indices.size() < NumEntries_[LocalRow]), std::runtime_error,
"Ifpack2::SingletonFilter::getLocalRowCopy invalid row or array size.");
309 LocalOrdinal ARow = InvReorder_[LocalRow];
310 A_->getLocalRowCopy(ARow,Indices_,Values_,Nnz);
314 for (
size_t i = 0 ; i < Nnz ; ++i) {
315 LocalOrdinal ii = Reorder_[Indices_[i]];
317 Indices[NumEntries] = ii;
318 Values[NumEntries] = Values_[i];
325template<
class MatrixType>
327 global_inds_host_view_type &,
328 values_host_view_type &)
const
330 throw std::runtime_error(
"Ifpack2::SingletonFilter: does not support getGlobalRowView.");
334template<
class MatrixType>
336 local_inds_host_view_type & ,
337 values_host_view_type & )
const
339 throw std::runtime_error(
"Ifpack2::SingletonFilter: does not support getLocalRowView.");
343template<
class MatrixType>
347 return A_->getLocalDiagCopy(diag);
350template<
class MatrixType>
353 throw std::runtime_error(
"Ifpack2::SingletonFilter does not support leftScale.");
356template<
class MatrixType>
359 throw std::runtime_error(
"Ifpack2::SingletonFilter does not support rightScale.");
362template<
class MatrixType>
364 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y,
365 Teuchos::ETransp mode,
369 typedef Scalar DomainScalar;
370 typedef Scalar RangeScalar;
374 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
375 "Ifpack2::SingletonFilter::apply ERROR: X.getNumVectors() != Y.getNumVectors().");
377 RangeScalar zero = Teuchos::ScalarTraits<RangeScalar>::zero();
378 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > x_ptr = X.get2dView();
379 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > y_ptr = Y.get2dViewNonConst();
382 size_t NumVectors = Y.getNumVectors();
385 for (
size_t i = 0 ; i < NumRows_ ; ++i) {
388 getLocalRowCopy(i,Indices_,Values_,Nnz);
389 if (mode==Teuchos::NO_TRANS){
390 for (
size_t j = 0 ; j < Nnz ; ++j)
391 for (
size_t k = 0 ; k < NumVectors ; ++k)
392 y_ptr[k][i] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][Indices_[j]];
394 else if (mode==Teuchos::TRANS){
395 for (
size_t j = 0 ; j < Nnz ; ++j)
396 for (
size_t k = 0 ; k < NumVectors ; ++k)
397 y_ptr[k][Indices_[j]] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][i];
400 for (
size_t j = 0 ; j < Nnz ; ++j)
401 for (
size_t k = 0 ; k < NumVectors ; ++k)
402 y_ptr[k][Indices_[j]] += Teuchos::ScalarTraits<RangeScalar>::conjugate((RangeScalar)Values_[j]) * (RangeScalar)x_ptr[k][i];
407template<
class MatrixType>
413template<
class MatrixType>
419template<
class MatrixType>
421 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
423 this->
template SolveSingletonsTempl<Scalar,Scalar>(RHS, LHS);
426template<
class MatrixType>
427template<
class DomainScalar,
class RangeScalar>
429 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
431 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > RHS_ptr = RHS.get2dView();
432 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > LHS_ptr = LHS.get2dViewNonConst();
434 for (
size_t i = 0 ; i < NumSingletons_ ; ++i) {
435 LocalOrdinal ii = SingletonIndex_[i];
439 for (
size_t j = 0 ; j < Nnz ; ++j) {
440 if (Indices_[j] == ii) {
441 for (
size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
442 LHS_ptr[k][ii] = (RangeScalar)RHS_ptr[k][ii] / (RangeScalar)Values_[j];
448template<
class MatrixType>
450 const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
451 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
453 this->
template CreateReducedRHSTempl<Scalar,Scalar>(LHS, RHS, ReducedRHS);
456template<
class MatrixType>
457template<
class DomainScalar,
class RangeScalar>
459 const Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
460 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
462 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const RangeScalar > > RHS_ptr = RHS.get2dView();
463 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > LHS_ptr = LHS.get2dView();
464 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > ReducedRHS_ptr = ReducedRHS.get2dViewNonConst();
466 size_t NumVectors = LHS.getNumVectors();
468 for (
size_t i = 0 ; i < NumRows_ ; ++i)
469 for (
size_t k = 0 ; k < NumVectors ; ++k)
470 ReducedRHS_ptr[k][i] = RHS_ptr[k][InvReorder_[i]];
472 for (
size_t i = 0 ; i < NumRows_ ; ++i) {
473 LocalOrdinal ii = InvReorder_[i];
477 for (
size_t j = 0 ; j < Nnz ; ++j) {
478 if (Reorder_[Indices_[j]] == -1) {
479 for (
size_t k = 0 ; k < NumVectors ; ++k)
480 ReducedRHS_ptr[k][i] -= (RangeScalar)Values_[j] * (RangeScalar)LHS_ptr[k][Indices_[j]];
486template<
class MatrixType>
488 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
490 this->
template UpdateLHSTempl<Scalar,Scalar>(ReducedLHS, LHS);
493template<
class MatrixType>
494template<
class DomainScalar,
class RangeScalar>
496 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
499 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > LHS_ptr = LHS.get2dViewNonConst();
500 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > ReducedLHS_ptr = ReducedLHS.get2dView();
502 for (
size_t i = 0 ; i < NumRows_ ; ++i)
503 for (
size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
504 LHS_ptr[k][InvReorder_[i]] = (RangeScalar)ReducedLHS_ptr[k][i];
507template<
class MatrixType>
510 throw std::runtime_error(
"Ifpack2::SingletonFilter does not implement getFrobeniusNorm.");
515#define IFPACK2_SINGLETONFILTER_INSTANT(S,LO,GO,N) \
516 template class Ifpack2::SingletonFilter< Tpetra::RowMatrix<S, LO, GO, N> >;
Filter based on matrix entries.
Definition: Ifpack2_SingletonFilter_decl.hpp:65
SingletonFilter(const Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Matrix)
Constructor.
Definition: Ifpack2_SingletonFilter_def.hpp:56
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:351
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 size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:235
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:420
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:344
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:335
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 size_t getLocalNumCols() const
Returns the number of columns needed to apply the forward operator on this node, i....
Definition: Ifpack2_SingletonFilter_def.hpp:217
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:277
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 global_size_t getGlobalNumCols() const
Returns the number of global columns in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:205
virtual mag_type getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:508
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
virtual global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:229
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:357
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:290
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 hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Definition: Ifpack2_SingletonFilter_def.hpp:408
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:301
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 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 bool hasColMap() const
Indicates whether this matrix has a well-defined column map.
Definition: Ifpack2_SingletonFilter_def.hpp:265
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:487
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:449
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:271
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 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:363
virtual global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:199
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 GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:223
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:326
virtual bool isFillComplete() const
Returns true if fillComplete() has been called.
Definition: Ifpack2_SingletonFilter_def.hpp:283
virtual ~SingletonFilter()
Destructor.
Definition: Ifpack2_SingletonFilter_def.hpp:143
virtual size_t getLocalNumRows() const
Returns the number of rows owned on the calling node.
Definition: Ifpack2_SingletonFilter_def.hpp:211
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the communicator.
Definition: Ifpack2_SingletonFilter_def.hpp:147
virtual bool supportsRowViews() const
Returns true if RowViews are supported.
Definition: Ifpack2_SingletonFilter_def.hpp:414
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:74