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" 56 template<
class MatrixType>
67 if (A_->getComm()->getSize() != 1 || A_->getNodeNumRows() != A_->getGlobalNumRows()) {
68 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.");
72 size_t NumRowsA_ = A_->getNodeNumRows();
76 MaxNumEntriesA_ = A_->getNodeMaxNumRowEntries();
79 Indices_.
resize(MaxNumEntriesA_);
80 Values_.
resize(MaxNumEntriesA_);
83 Reorder_.resize(NumRowsA_);
84 Reorder_.assign(Reorder_.size(),-1);
88 for (
size_t i = 0 ; i < NumRowsA_ ; ++i) {
90 A_->getLocalRowCopy(i,Indices_,Values_,Nnz);
92 Reorder_[i] = NumRows_++;
100 InvReorder_.resize(NumRows_);
101 for (
size_t i = 0 ; i < NumRowsA_ ; ++i) {
104 InvReorder_[Reorder_[i]] = i;
106 NumEntries_.resize(NumRows_);
107 SingletonIndex_.resize(NumSingletons_);
112 for (
size_t i = 0 ; i < NumRowsA_ ; ++i) {
114 A_->getLocalRowCopy(i,Indices_,Values_,Nnz);
115 LocalOrdinal ii = Reorder_[i];
117 NumEntries_[ii] = Nnz;
119 if (Nnz > MaxNumEntries_)
120 MaxNumEntries_ = Nnz;
123 SingletonIndex_[count] = i;
129 ReducedMap_ = Teuchos::rcp(
new Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node>(NumRows_,0,A_->getComm()) );
132 Diagonal_ = Teuchos::rcp(
new Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node>(ReducedMap_) );
134 Tpetra::Vector<Scalar,LocalOrdinal,GlobalOrdinal,Node> DiagonalA(A_->getRowMap());
135 A_->getLocalDiagCopy(DiagonalA);
137 for (
size_t i = 0 ; i < NumRows_ ; ++i) {
138 LocalOrdinal ii = InvReorder_[i];
139 Diagonal_->replaceLocalValue((LocalOrdinal)i,DiagonalAview[ii]);
144 template<
class MatrixType>
148 template<
class MatrixType>
152 return A_->getComm();
156 template<
class MatrixType>
160 return A_->getNode();
164 template<
class MatrixType>
165 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
166 typename MatrixType::global_ordinal_type,
167 typename MatrixType::node_type> >
174 template<
class MatrixType>
175 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
176 typename MatrixType::global_ordinal_type,
177 typename MatrixType::node_type> >
184 template<
class MatrixType>
185 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
186 typename MatrixType::global_ordinal_type,
187 typename MatrixType::node_type> >
194 template<
class MatrixType>
195 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
196 typename MatrixType::global_ordinal_type,
197 typename MatrixType::node_type> >
204 template<
class MatrixType>
205 Teuchos::RCP<
const Tpetra::RowGraph<
typename MatrixType::local_ordinal_type,
206 typename MatrixType::global_ordinal_type,
207 typename MatrixType::node_type> >
210 throw std::runtime_error(
"Ifpack2::SingletonFilter: does not support getGraph.");
214 template<
class MatrixType>
221 template<
class MatrixType>
228 template<
class MatrixType>
236 template<
class MatrixType>
243 template<
class MatrixType>
246 return A_->getIndexBase();
250 template<
class MatrixType>
257 template<
class MatrixType>
264 template<
class MatrixType>
267 throw std::runtime_error(
"Ifpack2::SingletonFilter does not implement getNumEntriesInGlobalRow.");
271 template<
class MatrixType>
274 return NumEntries_[localRow];
278 template<
class MatrixType>
281 return A_->getGlobalNumDiags();
285 template<
class MatrixType>
288 return A_->getNodeNumDiags();
292 template<
class MatrixType>
295 return MaxNumEntries_;
299 template<
class MatrixType>
302 return MaxNumEntries_;
306 template<
class MatrixType>
313 template<
class MatrixType>
316 return A_->isLowerTriangular();
320 template<
class MatrixType>
323 return A_->isUpperTriangular();
327 template<
class MatrixType>
330 return A_->isLocallyIndexed();
334 template<
class MatrixType>
337 return A_->isGloballyIndexed();
341 template<
class MatrixType>
344 return A_->isFillComplete();
348 template<
class MatrixType>
352 size_t &NumEntries)
const 354 throw std::runtime_error(
"Ifpack2::SingletonFilter does not implement getGlobalRowCopy.");
358 template<
class MatrixType>
362 size_t &NumEntries)
const 364 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.");
367 LocalOrdinal ARow = InvReorder_[LocalRow];
368 A_->getLocalRowCopy(ARow,Indices_(),Values_(),Nnz);
372 for (
size_t i = 0 ; i < Nnz ; ++i) {
373 LocalOrdinal ii = Reorder_[Indices_[i]];
375 Indices[NumEntries] = ii;
376 Values[NumEntries] = Values_[i];
384 template<
class MatrixType>
389 throw std::runtime_error(
"Ifpack2::SingletonFilter: does not support getGlobalRowView.");
393 template<
class MatrixType>
398 throw std::runtime_error(
"Ifpack2::SingletonFilter: does not support getLocalRowView.");
402 template<
class MatrixType>
406 return A_->getLocalDiagCopy(diag);
410 template<
class MatrixType>
413 throw std::runtime_error(
"Ifpack2::SingletonFilter does not support leftScale.");
417 template<
class MatrixType>
420 throw std::runtime_error(
"Ifpack2::SingletonFilter does not support rightScale.");
424 template<
class MatrixType>
426 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y,
431 typedef Scalar DomainScalar;
432 typedef Scalar RangeScalar;
437 "Ifpack2::SingletonFilter::apply ERROR: X.getNumVectors() != Y.getNumVectors().");
444 size_t NumVectors = Y.getNumVectors();
447 for (
size_t i = 0 ; i < NumRows_ ; ++i) {
450 getLocalRowCopy(i,Indices_(),Values_(),Nnz);
451 if (mode==Teuchos::NO_TRANS){
452 for (
size_t j = 0 ; j < Nnz ; ++j)
453 for (
size_t k = 0 ; k < NumVectors ; ++k)
454 y_ptr[k][i] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][Indices_[j]];
456 else if (mode==Teuchos::TRANS){
457 for (
size_t j = 0 ; j < Nnz ; ++j)
458 for (
size_t k = 0 ; k < NumVectors ; ++k)
459 y_ptr[k][Indices_[j]] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][i];
462 for (
size_t j = 0 ; j < Nnz ; ++j)
463 for (
size_t k = 0 ; k < NumVectors ; ++k)
471 template<
class MatrixType>
478 template<
class MatrixType>
485 template<
class MatrixType>
487 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
489 this->
template SolveSingletonsTempl<Scalar,Scalar>(RHS, LHS);
493 template<
class MatrixType>
494 template<
class DomainScalar,
class RangeScalar>
496 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
501 for (
size_t i = 0 ; i < NumSingletons_ ; ++i) {
502 LocalOrdinal ii = SingletonIndex_[i];
505 A_->getLocalRowCopy(ii,Indices_(),Values_(),Nnz);
506 for (
size_t j = 0 ; j < Nnz ; ++j) {
507 if (Indices_[j] == ii) {
508 for (
size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
509 LHS_ptr[k][ii] = (RangeScalar)RHS_ptr[k][ii] / (RangeScalar)Values_[j];
516 template<
class MatrixType>
518 const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
519 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
521 this->
template CreateReducedRHSTempl<Scalar,Scalar>(LHS, RHS, ReducedRHS);
525 template<
class MatrixType>
526 template<
class DomainScalar,
class RangeScalar>
528 const Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
529 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
535 size_t NumVectors = LHS.getNumVectors();
537 for (
size_t i = 0 ; i < NumRows_ ; ++i)
538 for (
size_t k = 0 ; k < NumVectors ; ++k)
539 ReducedRHS_ptr[k][i] = RHS_ptr[k][InvReorder_[i]];
541 for (
size_t i = 0 ; i < NumRows_ ; ++i) {
542 LocalOrdinal ii = InvReorder_[i];
544 A_->getLocalRowCopy(ii,Indices_(),Values_(),Nnz);
546 for (
size_t j = 0 ; j < Nnz ; ++j) {
547 if (Reorder_[Indices_[j]] == -1) {
548 for (
size_t k = 0 ; k < NumVectors ; ++k)
549 ReducedRHS_ptr[k][i] -= (RangeScalar)Values_[j] * (RangeScalar)LHS_ptr[k][Indices_[j]];
556 template<
class MatrixType>
558 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
560 this->
template UpdateLHSTempl<Scalar,Scalar>(ReducedLHS, LHS);
564 template<
class MatrixType>
565 template<
class DomainScalar,
class RangeScalar>
567 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
573 for (
size_t i = 0 ; i < NumRows_ ; ++i)
574 for (
size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
575 LHS_ptr[k][InvReorder_[i]] = (RangeScalar)ReducedLHS_ptr[k][i];
579 template<
class MatrixType>
582 throw std::runtime_error(
"Ifpack2::SingletonFilter does not implement getFrobeniusNorm.");
587 #define IFPACK2_SINGLETONFILTER_INSTANT(S,LO,GO,N) \ 588 template class Ifpack2::SingletonFilter< Tpetra::RowMatrix<S, LO, GO, N> >; 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:418
virtual ~SingletonFilter()
Destructor.
Definition: Ifpack2_SingletonFilter_def.hpp:145
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:328
virtual global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:215
virtual size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
Definition: Ifpack2_SingletonFilter_def.hpp:293
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:198
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:557
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:517
virtual void getGlobalRowView(GlobalOrdinal GlobalRow, Teuchos::ArrayView< const GlobalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:385
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual size_t getNodeNumDiags() const
Returns the number of local diagonal entries, based on global row/column index comparisons.
Definition: Ifpack2_SingletonFilter_def.hpp:286
virtual void getLocalRowCopy(LocalOrdinal LocalRow, const Teuchos::ArrayView< LocalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &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:359
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:188
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:168
virtual size_t getNodeNumRows() const
Returns the number of rows owned on the calling node.
Definition: Ifpack2_SingletonFilter_def.hpp:229
virtual bool hasTransposeApply() const
Indicates whether this operator supports applying the adjoint operator.
Definition: Ifpack2_SingletonFilter_def.hpp:472
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the communicator.
Definition: Ifpack2_SingletonFilter_def.hpp:150
virtual mag_type getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:580
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:265
virtual size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:258
virtual global_size_t getGlobalNumCols() const
Returns the number of global columns in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:222
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:425
LinearOp zero(const VectorSpace &vs)
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:403
void resize(size_type new_size, const value_type &x=value_type())
virtual bool isLowerTriangular() const
Indicates whether this matrix is lower triangular.
Definition: Ifpack2_SingletonFilter_def.hpp:314
virtual bool hasColMap() const
Indicates whether this matrix has a well-defined column map.
Definition: Ifpack2_SingletonFilter_def.hpp:307
virtual Teuchos::RCP< const Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node > > getGraph() const
Returns the RowGraph associated with this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:208
virtual bool isUpperTriangular() const
Indicates whether this matrix is upper triangular.
Definition: Ifpack2_SingletonFilter_def.hpp:321
virtual size_t getNodeNumCols() 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:237
virtual global_size_t getGlobalNumDiags() const
Returns the number of global diagonal entries, based on global row/column index comparisons.
Definition: Ifpack2_SingletonFilter_def.hpp:279
SingletonFilter(const Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &Matrix)
Constructor.
Definition: Ifpack2_SingletonFilter_def.hpp:57
virtual Teuchos::RCP< Node > getNode() const
Returns the underlying node.
Definition: Ifpack2_SingletonFilter_def.hpp:158
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:272
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:486
virtual bool supportsRowViews() const
Returns true if RowViews are supported.
Definition: Ifpack2_SingletonFilter_def.hpp:479
virtual void getLocalRowView(LocalOrdinal LocalRow, Teuchos::ArrayView< const LocalOrdinal > &indices, Teuchos::ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:394
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:72
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:178
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:411
virtual bool isFillComplete() const
Returns true if fillComplete() has been called.
Definition: Ifpack2_SingletonFilter_def.hpp:342
virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow, const Teuchos::ArrayView< GlobalOrdinal > &Indices, const Teuchos::ArrayView< Scalar > &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:349
virtual size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
Definition: Ifpack2_SingletonFilter_def.hpp:300
virtual global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:251
virtual GlobalOrdinal getIndexBase() const
Returns the index base for global indices for this matrix.
Definition: Ifpack2_SingletonFilter_def.hpp:244
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:335