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"
55 template<
class MatrixType>
66 if (A_->getComm()->getSize() != 1 || A_->getNodeNumRows() != 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_->getNodeNumRows();
75 MaxNumEntriesA_ = A_->getNodeMaxNumRowEntries();
78 Indices_.resize(MaxNumEntriesA_);
79 Values_.resize(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]);
142 template<
class MatrixType>
145 template<
class MatrixType>
146 Teuchos::RCP<const Teuchos::Comm<int> >
149 return A_->getComm();
152 template<
class MatrixType>
153 Teuchos::RCP<typename MatrixType::node_type>
156 return A_->getNode();
159 template<
class MatrixType>
160 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
161 typename MatrixType::global_ordinal_type,
162 typename MatrixType::node_type> >
168 template<
class MatrixType>
169 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
170 typename MatrixType::global_ordinal_type,
171 typename MatrixType::node_type> >
177 template<
class MatrixType>
178 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
179 typename MatrixType::global_ordinal_type,
180 typename MatrixType::node_type> >
186 template<
class MatrixType>
187 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
188 typename MatrixType::global_ordinal_type,
189 typename MatrixType::node_type> >
195 template<
class MatrixType>
196 Teuchos::RCP<
const Tpetra::RowGraph<
typename MatrixType::local_ordinal_type,
197 typename MatrixType::global_ordinal_type,
198 typename MatrixType::node_type> >
201 throw std::runtime_error(
"Ifpack2::SingletonFilter: does not support getGraph.");
204 template<
class MatrixType>
210 template<
class MatrixType>
216 template<
class MatrixType>
222 template<
class MatrixType>
228 template<
class MatrixType>
231 return A_->getIndexBase();
234 template<
class MatrixType>
240 template<
class MatrixType>
246 template<
class MatrixType>
249 throw std::runtime_error(
"Ifpack2::SingletonFilter does not implement getNumEntriesInGlobalRow.");
252 template<
class MatrixType>
255 return NumEntries_[localRow];
258 template<
class MatrixType>
261 return MaxNumEntries_;
264 template<
class MatrixType>
267 return MaxNumEntries_;
270 template<
class MatrixType>
276 template<
class MatrixType>
279 return A_->isLocallyIndexed();
282 template<
class MatrixType>
285 return A_->isGloballyIndexed();
288 template<
class MatrixType>
291 return A_->isFillComplete();
294 template<
class MatrixType>
296 const Teuchos::ArrayView<GlobalOrdinal> &,
297 const Teuchos::ArrayView<Scalar> &,
300 throw std::runtime_error(
"Ifpack2::SingletonFilter does not implement getGlobalRowCopy.");
303 template<
class MatrixType>
305 const Teuchos::ArrayView<LocalOrdinal> &Indices,
306 const Teuchos::ArrayView<Scalar> &Values,
307 size_t &NumEntries)
const
309 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.");
312 LocalOrdinal ARow = InvReorder_[LocalRow];
313 A_->getLocalRowCopy(ARow,Indices_(),Values_(),Nnz);
317 for (
size_t i = 0 ; i < Nnz ; ++i) {
318 LocalOrdinal ii = Reorder_[Indices_[i]];
320 Indices[NumEntries] = ii;
321 Values[NumEntries] = Values_[i];
328 template<
class MatrixType>
330 Teuchos::ArrayView<const GlobalOrdinal> &,
331 Teuchos::ArrayView<const Scalar> &)
const
333 throw std::runtime_error(
"Ifpack2::SingletonFilter: does not support getGlobalRowView.");
336 template<
class MatrixType>
338 Teuchos::ArrayView<const LocalOrdinal> &,
339 Teuchos::ArrayView<const Scalar> &)
const
341 throw std::runtime_error(
"Ifpack2::SingletonFilter: does not support getLocalRowView.");
344 template<
class MatrixType>
348 return A_->getLocalDiagCopy(diag);
351 template<
class MatrixType>
354 throw std::runtime_error(
"Ifpack2::SingletonFilter does not support leftScale.");
357 template<
class MatrixType>
360 throw std::runtime_error(
"Ifpack2::SingletonFilter does not support rightScale.");
363 template<
class MatrixType>
365 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y,
366 Teuchos::ETransp mode,
370 typedef Scalar DomainScalar;
371 typedef Scalar RangeScalar;
375 TEUCHOS_TEST_FOR_EXCEPTION(X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
376 "Ifpack2::SingletonFilter::apply ERROR: X.getNumVectors() != Y.getNumVectors().");
378 RangeScalar zero = Teuchos::ScalarTraits<RangeScalar>::zero();
379 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > x_ptr = X.get2dView();
380 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > y_ptr = Y.get2dViewNonConst();
383 size_t NumVectors = Y.getNumVectors();
386 for (
size_t i = 0 ; i < NumRows_ ; ++i) {
389 getLocalRowCopy(i,Indices_(),Values_(),Nnz);
390 if (mode==Teuchos::NO_TRANS){
391 for (
size_t j = 0 ; j < Nnz ; ++j)
392 for (
size_t k = 0 ; k < NumVectors ; ++k)
393 y_ptr[k][i] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][Indices_[j]];
395 else if (mode==Teuchos::TRANS){
396 for (
size_t j = 0 ; j < Nnz ; ++j)
397 for (
size_t k = 0 ; k < NumVectors ; ++k)
398 y_ptr[k][Indices_[j]] += (RangeScalar)Values_[j] * (RangeScalar)x_ptr[k][i];
401 for (
size_t j = 0 ; j < Nnz ; ++j)
402 for (
size_t k = 0 ; k < NumVectors ; ++k)
403 y_ptr[k][Indices_[j]] += Teuchos::ScalarTraits<RangeScalar>::conjugate((RangeScalar)Values_[j]) * (RangeScalar)x_ptr[k][i];
408 template<
class MatrixType>
414 template<
class MatrixType>
420 template<
class MatrixType>
422 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
424 this->
template SolveSingletonsTempl<Scalar,Scalar>(RHS, LHS);
427 template<
class MatrixType>
428 template<
class DomainScalar,
class RangeScalar>
430 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
432 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > RHS_ptr = RHS.get2dView();
433 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > LHS_ptr = LHS.get2dViewNonConst();
435 for (
size_t i = 0 ; i < NumSingletons_ ; ++i) {
436 LocalOrdinal ii = SingletonIndex_[i];
440 for (
size_t j = 0 ; j < Nnz ; ++j) {
441 if (Indices_[j] == ii) {
442 for (
size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
443 LHS_ptr[k][ii] = (RangeScalar)RHS_ptr[k][ii] / (RangeScalar)Values_[j];
449 template<
class MatrixType>
451 const Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
452 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
454 this->
template CreateReducedRHSTempl<Scalar,Scalar>(LHS, RHS, ReducedRHS);
457 template<
class MatrixType>
458 template<
class DomainScalar,
class RangeScalar>
460 const Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& RHS,
461 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& ReducedRHS)
463 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const RangeScalar > > RHS_ptr = RHS.get2dView();
464 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > LHS_ptr = LHS.get2dView();
465 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > ReducedRHS_ptr = ReducedRHS.get2dViewNonConst();
467 size_t NumVectors = LHS.getNumVectors();
469 for (
size_t i = 0 ; i < NumRows_ ; ++i)
470 for (
size_t k = 0 ; k < NumVectors ; ++k)
471 ReducedRHS_ptr[k][i] = RHS_ptr[k][InvReorder_[i]];
473 for (
size_t i = 0 ; i < NumRows_ ; ++i) {
474 LocalOrdinal ii = InvReorder_[i];
478 for (
size_t j = 0 ; j < Nnz ; ++j) {
479 if (Reorder_[Indices_[j]] == -1) {
480 for (
size_t k = 0 ; k < NumVectors ; ++k)
481 ReducedRHS_ptr[k][i] -= (RangeScalar)Values_[j] * (RangeScalar)LHS_ptr[k][Indices_[j]];
487 template<
class MatrixType>
489 Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
491 this->
template UpdateLHSTempl<Scalar,Scalar>(ReducedLHS, LHS);
494 template<
class MatrixType>
495 template<
class DomainScalar,
class RangeScalar>
497 Tpetra::MultiVector<RangeScalar,LocalOrdinal,GlobalOrdinal,Node>& LHS)
500 Teuchos::ArrayRCP<Teuchos::ArrayRCP<RangeScalar> > LHS_ptr = LHS.get2dViewNonConst();
501 Teuchos::ArrayRCP<Teuchos::ArrayRCP<const DomainScalar> > ReducedLHS_ptr = ReducedLHS.get2dView();
503 for (
size_t i = 0 ; i < NumRows_ ; ++i)
504 for (
size_t k = 0 ; k < LHS.getNumVectors() ; ++k)
505 LHS_ptr[k][InvReorder_[i]] = (RangeScalar)ReducedLHS_ptr[k][i];
508 template<
class MatrixType>
511 throw std::runtime_error(
"Ifpack2::SingletonFilter does not implement getFrobeniusNorm.");
516 #define IFPACK2_SINGLETONFILTER_INSTANT(S,LO,GO,N) \
517 template class Ifpack2::SingletonFilter< Tpetra::RowMatrix<S, LO, GO, N> >;