43 #ifndef IFPACK2_REORDERFILTER_DEF_HPP
44 #define IFPACK2_REORDERFILTER_DEF_HPP
45 #include "Ifpack2_ReorderFilter_decl.hpp"
48 #include "Tpetra_ConfigDefs.hpp"
49 #include "Tpetra_RowMatrix.hpp"
50 #include "Tpetra_Map.hpp"
51 #include "Tpetra_MultiVector.hpp"
52 #include "Tpetra_Vector.hpp"
56 template<
class MatrixType>
63 reverseperm_ (reverseperm)
66 A_.
is_null (), std::invalid_argument,
67 "Ifpack2::ReorderFilter: The input matrix is null.");
71 A_->getComm()->getSize() != 1, std::invalid_argument,
72 "Ifpack2::ReorderFilter: This class may only be used if the input matrix's "
73 "communicator has one process. This class is an implementation detail of "
74 "Ifpack2::AdditiveSchwarz, and it is not meant to be used otherwise.");
77 A_->getLocalNumRows () != A_->getGlobalNumRows (),
78 std::invalid_argument,
79 "Ifpack2::ReorderFilter: The input matrix is not square.");
82 Kokkos::resize(Indices_,A_->getLocalMaxNumRowEntries ());
83 Kokkos::resize(Values_,A_->getLocalMaxNumRowEntries ());
87 template<
class MatrixType>
91 template<
class MatrixType>
100 template<
class MatrixType>
105 A_.
is_null (), std::runtime_error,
"Ifpack2::ReorderFilter::"
106 "getRowMap: The matrix A is null, so there is no row Map.");
108 return A_->getRowMap ();
112 template<
class MatrixType>
117 A_.
is_null (), std::runtime_error,
"Ifpack2::ReorderFilter::"
118 "getColMap: The matrix A is null, so there is no column Map.");
120 return A_->getColMap();
124 template<
class MatrixType>
129 A_.
is_null (), std::runtime_error,
"Ifpack2::ReorderFilter::"
130 "getDomainMap: The matrix A is null, so there is no domain Map.");
132 return A_->getDomainMap();
136 template<
class MatrixType>
141 A_.
is_null (), std::runtime_error,
"Ifpack2::ReorderFilter::"
142 "getRangeMap: The matrix A is null, so there is no range Map.");
144 return A_->getRangeMap();
148 template<
class MatrixType>
149 Teuchos::RCP<
const Tpetra::RowGraph<
typename MatrixType::local_ordinal_type,
150 typename MatrixType::global_ordinal_type,
151 typename MatrixType::node_type> >
154 throw std::runtime_error(
"Ifpack2::ReorderFilter: does not support getGraph.");
158 template<
class MatrixType>
161 return A_->getGlobalNumRows();
165 template<
class MatrixType>
168 return A_->getGlobalNumCols();
172 template<
class MatrixType>
175 return A_->getLocalNumRows();
179 template<
class MatrixType>
182 return A_->getLocalNumCols();
186 template<
class MatrixType>
189 return A_->getIndexBase();
193 template<
class MatrixType>
196 return A_->getGlobalNumEntries();
200 template<
class MatrixType>
203 return A_->getLocalNumEntries();
206 template<
class MatrixType>
209 return A_->getBlockSize();
212 template<
class MatrixType>
220 const local_ordinal_type lclRow =
221 A_->getRowMap ()->getLocalElement (globalRow);
224 return static_cast<size_t> (0);
226 const local_ordinal_type origLclRow = reverseperm_[lclRow];
227 return A_->getNumEntriesInLocalRow (origLclRow);
232 template<
class MatrixType>
238 if (A_->getRowMap ()->isNodeLocalElement (localRow)) {
240 const local_ordinal_type localReorderedRow = reverseperm_[localRow];
241 return A_->getNumEntriesInLocalRow (localReorderedRow);
244 return static_cast<size_t> (0);
249 template<
class MatrixType>
252 return A_->getGlobalMaxNumRowEntries();
256 template<
class MatrixType>
259 return A_->getLocalMaxNumRowEntries();
263 template<
class MatrixType>
270 template<
class MatrixType>
273 return A_->isLocallyIndexed();
277 template<
class MatrixType>
280 return A_->isGloballyIndexed();
284 template<
class MatrixType>
287 return A_->isFillComplete();
291 template<
class MatrixType>
294 nonconst_global_inds_host_view_type &globalInd,
295 nonconst_values_host_view_type &val,
296 size_t& numEntries)
const
300 using Teuchos::av_reinterpret_cast;
301 typedef local_ordinal_type LO;
304 const map_type& rowMap = * (A_->getRowMap ());
305 const local_ordinal_type localRow = rowMap.getLocalElement (globalRow);
307 localRow == OTLO::invalid (), std::invalid_argument,
"Ifpack2::Reorder"
308 "Filter::getGlobalRowCopy: The given global row index " << globalRow
309 <<
" is not owned by the calling process with rank "
310 << rowMap.getComm ()->getRank () <<
".");
313 numEntries = this->getNumEntriesInLocalRow (localRow);
314 this->getLocalRowCopy (localRow, Indices_, val, numEntries);
317 for (
size_t k = 0; k < numEntries; ++k) {
318 globalInd[k] = rowMap.getGlobalElement (Indices_[k]);
323 template<
class MatrixType>
326 nonconst_local_inds_host_view_type &Indices,
327 nonconst_values_host_view_type &Values,
328 size_t& NumEntries)
const
332 ! A_->getRowMap ()->isNodeLocalElement (LocalRow),
333 std::invalid_argument,
334 "Ifpack2::ReorderFilter::getLocalRowCopy: The given local row index "
335 << LocalRow <<
" is not a valid local row index on the calling process "
336 "with rank " << A_->getRowMap ()->getComm ()->getRank () <<
".");
340 const local_ordinal_type origLclRow = reverseperm_[LocalRow];
341 const size_t numEntries = A_->getNumEntriesInLocalRow (origLclRow);
344 static_cast<size_t> (Indices.size ()) < numEntries ||
345 static_cast<size_t> (Values.size ()) < numEntries,
346 std::invalid_argument,
347 "Ifpack2::ReorderFilter::getLocalRowCopy: The given array views are not "
348 "long enough to store all the data in the given row " << LocalRow
349 <<
". Indices.size() = " << Indices.size () <<
", Values.size() = "
350 << Values.size () <<
", but the (original) row has " << numEntries
353 A_->getLocalRowCopy (origLclRow, Indices, Values, NumEntries);
358 for (
size_t i = 0; i < NumEntries; ++i) {
359 Indices[i] = perm_[Indices[i]];
364 template<
class MatrixType>
366 global_inds_host_view_type &,
367 values_host_view_type &)
const
369 throw std::runtime_error(
"Ifpack2::ReorderFilter: does not support getGlobalRowView.");
374 template<
class MatrixType>
376 local_inds_host_view_type & ,
377 values_host_view_type & )
const
379 throw std::runtime_error(
"Ifpack2::ReorderFilter: does not support getLocalRowView.");
384 template<
class MatrixType>
386 getLocalDiagCopy (Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &diag)
const
389 return A_->getLocalDiagCopy(diag);
393 template<
class MatrixType>
396 throw std::runtime_error(
"Ifpack2::ReorderFilter does not support leftScale.");
400 template<
class MatrixType>
403 throw std::runtime_error(
"Ifpack2::ReorderFilter does not support rightScale.");
407 template<
class MatrixType>
409 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &X,
410 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &Y,
413 scalar_type beta)
const
418 alpha != STS::one () || beta != STS::zero (), std::logic_error,
419 "Ifpack2::ReorderFilter::apply is only implemented for alpha = 1 and "
420 "beta = 0. You set alpha = " << alpha <<
" and beta = " << beta <<
".");
425 X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
426 "Ifpack2::ReorderFilter::apply: X.getNumVectors() != Y.getNumVectors().");
428 const scalar_type zero = STS::zero ();
433 const size_t NumVectors = Y.getNumVectors ();
435 for (
size_t i = 0; i < A_->getLocalNumRows (); ++i) {
438 getLocalRowCopy (i, Indices_ , Values_ , Nnz);
439 scalar_type* Values =
reinterpret_cast<scalar_type*
>(Values_.data());
441 for (
size_t j = 0; j < Nnz; ++j) {
442 for (
size_t k = 0; k < NumVectors; ++k) {
443 y_ptr[k][i] += Values[j] * x_ptr[k][Indices_[j]];
448 for (
size_t j = 0; j < Nnz; ++j) {
449 for (
size_t k = 0; k < NumVectors; ++k) {
450 y_ptr[k][Indices_[j]] += Values[j] * x_ptr[k][i];
455 for (
size_t j = 0; j < Nnz; ++j) {
456 for (
size_t k = 0; k < NumVectors; ++k) {
457 y_ptr[k][Indices_[j]] += STS::conjugate(Values[j]) * x_ptr[k][i];
465 template<
class MatrixType>
472 template<
class MatrixType>
479 template<
class MatrixType>
483 return A_->getFrobeniusNorm ();
487 template<
class MatrixType>
490 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &reorderedY)
const
492 this->
template permuteOriginalToReorderedTempl<scalar_type,scalar_type>(originalX, reorderedY);
496 template<
class MatrixType>
497 template<
class DomainScalar,
class RangeScalar>
499 Tpetra::MultiVector<RangeScalar,local_ordinal_type,global_ordinal_type,node_type> &reorderedY)
const
502 "Ifpack2::ReorderFilter::permuteOriginalToReordered ERROR: X.getNumVectors() != Y.getNumVectors().");
507 const local_ordinal_type blockSize = getBlockSize();
508 const local_ordinal_type numRows = originalX.getLocalLength() / blockSize;
509 for(
size_t k=0; k < originalX.getNumVectors(); k++)
510 for(local_ordinal_type i=0; i< numRows; i++)
511 for(local_ordinal_type j=0; j< blockSize; ++j)
512 y_ptr[k][perm_[i]*blockSize + j] = (RangeScalar)x_ptr[k][i*blockSize + j];
516 template<
class MatrixType>
518 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &originalY)
const
520 this->
template permuteReorderedToOriginalTempl<scalar_type,scalar_type>(reorderedX, originalY);
524 template<
class MatrixType>
525 template<
class DomainScalar,
class RangeScalar>
528 Tpetra::MultiVector<RangeScalar,local_ordinal_type,global_ordinal_type,node_type> &originalY)
const
531 reorderedX.getNumVectors() != originalY.getNumVectors(),
533 "Ifpack2::ReorderFilter::permuteReorderedToOriginal: "
534 "X.getNumVectors() != Y.getNumVectors().");
536 #ifdef HAVE_IFPACK2_DEBUG
540 reorderedX.norm2 (norms ());
543 j < reorderedX.getNumVectors (); ++j) {
544 if (STM::isnaninf (norms[j])) {
550 ! good, std::runtime_error,
"Ifpack2::ReorderFilter::"
551 "permuteReorderedToOriginalTempl: The 2-norm of the input reorderedX is "
554 #endif // HAVE_IFPACK2_DEBUG
559 const local_ordinal_type blockSize = getBlockSize();
560 const local_ordinal_type numRows = reorderedX.getLocalLength() / blockSize;
561 for (
size_t k = 0; k < reorderedX.getNumVectors (); ++k) {
562 for (local_ordinal_type i = 0; i < numRows; ++i) {
563 for(local_ordinal_type j = 0; j < blockSize; ++j) {
564 y_ptr[k][reverseperm_[i]*blockSize + j] = (RangeScalar) x_ptr[k][i*blockSize + j];
569 #ifdef HAVE_IFPACK2_DEBUG
573 originalY.norm2 (norms ());
576 j < originalY.getNumVectors (); ++j) {
577 if (STM::isnaninf (norms[j])) {
583 ! good, std::runtime_error,
"Ifpack2::ReorderFilter::"
584 "permuteReorderedToOriginalTempl: The 2-norm of the output originalY is "
587 #endif // HAVE_IFPACK2_DEBUG
592 #define IFPACK2_REORDERFILTER_INSTANT(S,LO,GO,N) \
593 template class Ifpack2::ReorderFilter< Tpetra::RowMatrix<S, LO, GO, N> >;
virtual size_t getLocalNumRows() const
Returns the number of rows owned on the calling node.
Definition: Ifpack2_ReorderFilter_def.hpp:173
Wraps a Tpetra::RowMatrix in a filter that reorders local rows and columns.
Definition: Ifpack2_ReorderFilter_decl.hpp:69
virtual mag_type getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:480
virtual Teuchos::RCP< const map_type > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:102
virtual bool supportsRowViews() const
Returns true if RowViews are supported.
Definition: Ifpack2_ReorderFilter_def.hpp:473
virtual size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
Definition: Ifpack2_ReorderFilter_def.hpp:250
virtual Teuchos::RCP< const map_type > getDomainMap() const
Returns the Map that describes the domain distribution in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:126
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:159
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_ReorderFilter_def.hpp:180
virtual Teuchos::RCP< const map_type > getColMap() const
Returns the Map that describes the column distribution in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:114
ReorderFilter(const Teuchos::RCP< const row_matrix_type > &A, const Teuchos::ArrayRCP< local_ordinal_type > &perm, const Teuchos::ArrayRCP< local_ordinal_type > &reverseperm)
Constructor.
Definition: Ifpack2_ReorderFilter_def.hpp:58
virtual Teuchos::RCP< const map_type > getRangeMap() const
Returns the Map that describes the range distribution in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:138
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_ReorderFilter_def.hpp:365
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
, where Op(A) is either A, , or .
Definition: Ifpack2_ReorderFilter_def.hpp:409
virtual void getGlobalRowCopy(global_ordinal_type 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_ReorderFilter_def.hpp:293
virtual bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Definition: Ifpack2_ReorderFilter_def.hpp:278
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The matrix's communicator.
Definition: Ifpack2_ReorderFilter_def.hpp:92
virtual void getLocalDiagCopy(Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
Definition: Ifpack2_ReorderFilter_def.hpp:386
virtual size_t getNumEntriesInLocalRow(local_ordinal_type localRow) const
The current number of entries in this matrix, stored on the calling process, in the row whose local i...
Definition: Ifpack2_ReorderFilter_def.hpp:234
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_ReorderFilter_def.hpp:375
virtual bool hasColMap() const
Indicates whether this matrix has a well-defined column map.
Definition: Ifpack2_ReorderFilter_def.hpp:264
virtual size_t getNumEntriesInGlobalRow(global_ordinal_type globalRow) const
The current number of entries in this matrix, stored on the calling process, in the row whose global ...
Definition: Ifpack2_ReorderFilter_def.hpp:214
virtual global_size_t getGlobalNumCols() const
Returns the number of global columns in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:166
virtual bool hasTransposeApply() const
Whether apply() can apply the transpose or conjugate transpose.
Definition: Ifpack2_ReorderFilter_def.hpp:466
virtual size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:201
virtual global_ordinal_type getIndexBase() const
Returns the index base for global indices for this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:187
virtual global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:194
virtual size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
Definition: Ifpack2_ReorderFilter_def.hpp:257
virtual bool isFillComplete() const
Returns true if fillComplete() has been called.
Definition: Ifpack2_ReorderFilter_def.hpp:285
virtual void permuteReorderedToOriginal(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &reorderedX, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &originalY) const
Permute multivector: reordered-to-original.
Definition: Ifpack2_ReorderFilter_def.hpp:517
virtual Teuchos::RCP< const Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > > getGraph() const
Returns the RowGraph associated with this matrix.
Definition: Ifpack2_ReorderFilter_def.hpp:152
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_ReorderFilter_def.hpp:394
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_ReorderFilter_def.hpp:401
virtual void getLocalRowCopy(local_ordinal_type 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_ReorderFilter_def.hpp:325
virtual void permuteOriginalToReordered(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &originalX, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &reorderedY) const
Permute multivector: original-to-reordered.
Definition: Ifpack2_ReorderFilter_def.hpp:489
virtual ~ReorderFilter()
Destructor.
Definition: Ifpack2_ReorderFilter_def.hpp:88
virtual bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
Definition: Ifpack2_ReorderFilter_def.hpp:271
virtual local_ordinal_type getBlockSize() const
The number of degrees of freedom per mesh point.
Definition: Ifpack2_ReorderFilter_def.hpp:207