43 #ifndef IFPACK2_LOCALFILTER_DEF_HPP
44 #define IFPACK2_LOCALFILTER_DEF_HPP
46 #include <Ifpack2_LocalFilter_decl.hpp>
47 #include <Tpetra_Map.hpp>
48 #include <Tpetra_MultiVector.hpp>
49 #include <Tpetra_Vector.hpp>
54 # include "Teuchos_DefaultSerialComm.hpp"
60 template<
class MatrixType>
62 LocalFilter<MatrixType>::
63 mapPairsAreFitted (
const row_matrix_type& A)
65 const map_type& rangeMap = * (A.getRangeMap ());
66 const map_type& rowMap = * (A.getRowMap ());
67 const bool rangeAndRowFitted = mapPairIsFitted (rangeMap, rowMap);
69 const map_type& domainMap = * (A.getDomainMap ());
70 const map_type& columnMap = * (A.getColMap ());
71 const bool domainAndColumnFitted = mapPairIsFitted (domainMap, columnMap);
73 return rangeAndRowFitted && domainAndColumnFitted;
77 template<
class MatrixType>
79 LocalFilter<MatrixType>::
80 mapPairIsFitted (
const map_type& map1,
const map_type& map2)
82 return map1.isLocallyFitted (map2);
86 template<
class MatrixType>
97 #ifdef HAVE_IFPACK2_DEBUG
99 ! mapPairsAreFitted (*A), std::invalid_argument,
"Ifpack2::LocalFilter: "
100 "A's Map pairs are not fitted to each other on Process "
101 << A_->getRowMap ()->getComm ()->getRank () <<
" of the input matrix's "
103 "This means that LocalFilter does not currently know how to work with A. "
104 "This will change soon. Please see discussion of Bug 5992.");
105 #endif // HAVE_IFPACK2_DEBUG
108 RCP<const Teuchos::Comm<int> > localComm;
110 localComm =
rcp (
new Teuchos::MpiComm<int> (MPI_COMM_SELF));
138 const size_t numRows = A_->getRangeMap()->getNodeNumElements ();
146 const global_ordinal_type indexBase =
static_cast<global_ordinal_type
> (0);
150 Tpetra::GloballyDistributed));
153 localRangeMap_ = localRowMap_;
160 localDomainMap_ = localRangeMap_;
163 const size_t numCols = A_->getDomainMap()->getNodeNumElements ();
166 Tpetra::GloballyDistributed));
172 NumEntries_.resize (numRows);
176 MaxNumEntries_ = A_->getNodeMaxNumRowEntries ();
177 MaxNumEntriesA_ = A_->getNodeMaxNumRowEntries ();
180 localIndices_.
resize (MaxNumEntries_);
181 Values_.
resize (MaxNumEntries_);
190 size_t ActualMaxNumEntries = 0;
192 for (
size_t i = 0; i < numRows; ++i) {
194 size_t Nnz, NewNnz = 0;
195 A_->getLocalRowCopy (i, localIndices_, Values_, Nnz);
196 for (
size_t j = 0; j < Nnz; ++j) {
208 if (static_cast<size_t> (localIndices_[j]) < numRows) {
213 if (NewNnz > ActualMaxNumEntries) {
214 ActualMaxNumEntries = NewNnz;
217 NumNonzeros_ += NewNnz;
218 NumEntries_[i] = NewNnz;
221 MaxNumEntries_ = ActualMaxNumEntries;
225 template<
class MatrixType>
230 template<
class MatrixType>
234 return localRowMap_->getComm ();
240 template<
class MatrixType>
241 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
242 typename MatrixType::global_ordinal_type,
243 typename MatrixType::node_type> >
250 template<
class MatrixType>
251 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
252 typename MatrixType::global_ordinal_type,
253 typename MatrixType::node_type> >
260 template<
class MatrixType>
261 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
262 typename MatrixType::global_ordinal_type,
263 typename MatrixType::node_type> >
266 return localDomainMap_;
270 template<
class MatrixType>
271 Teuchos::RCP<
const Tpetra::Map<
typename MatrixType::local_ordinal_type,
272 typename MatrixType::global_ordinal_type,
273 typename MatrixType::node_type> >
276 return localRangeMap_;
280 template<
class MatrixType>
281 Teuchos::RCP<
const Tpetra::RowGraph<
typename MatrixType::local_ordinal_type,
282 typename MatrixType::global_ordinal_type,
283 typename MatrixType::node_type> >
289 return A_->getGraph ();
293 template<
class MatrixType>
296 return static_cast<global_size_t
> (localRangeMap_->getNodeNumElements ());
300 template<
class MatrixType>
303 return static_cast<global_size_t
> (localDomainMap_->getNodeNumElements ());
307 template<
class MatrixType>
310 return static_cast<size_t> (localRangeMap_->getNodeNumElements ());
314 template<
class MatrixType>
317 return static_cast<size_t> (localDomainMap_->getNodeNumElements ());
321 template<
class MatrixType>
322 typename MatrixType::global_ordinal_type
325 return A_->getIndexBase ();
329 template<
class MatrixType>
336 template<
class MatrixType>
343 template<
class MatrixType>
357 return NumEntries_[localRow];
362 template<
class MatrixType>
372 if (getRowMap ()->isNodeLocalElement (localRow)) {
373 return NumEntries_[localRow];
385 template<
class MatrixType>
388 return MaxNumEntries_;
392 template<
class MatrixType>
395 return MaxNumEntries_;
399 template<
class MatrixType>
406 template<
class MatrixType>
409 return A_->isLocallyIndexed ();
413 template<
class MatrixType>
416 return A_->isGloballyIndexed();
420 template<
class MatrixType>
423 return A_->isFillComplete ();
427 template<
class MatrixType>
433 size_t& numEntries)
const
438 const LO localRow = this->getRowMap ()->getLocalElement (globalRow);
451 numEntries = this->getNumEntriesInLocalRow (localRow);
456 this->getLocalRowCopy (localRow, localIndices (), values, numEntries);
458 const map_type& colMap = * (this->getColMap ());
461 const size_type numEnt =
462 std::min (static_cast<size_type> (numEntries),
463 std::min (globalIndices.
size (), values.
size ()));
464 for (size_type k = 0; k < numEnt; ++k) {
465 globalIndices[k] = colMap.getGlobalElement (localIndices[k]);
471 template<
class MatrixType>
477 size_t &NumEntries)
const
480 typedef global_ordinal_type GO;
482 if (! A_->getRowMap ()->isNodeLocalElement (LocalRow)) {
488 if (A_->getRowMap()->getComm()->getSize() == 1) {
489 A_->getLocalRowCopy (LocalRow, Indices (), Values (), NumEntries);
494 const size_t numEntInLclRow = NumEntries_[LocalRow];
495 if (static_cast<size_t> (Indices.
size ()) < numEntInLclRow ||
496 static_cast<size_t> (Values.
size ()) < numEntInLclRow) {
501 true, std::runtime_error,
502 "Ifpack2::LocalFilter::getLocalRowCopy: Invalid output array length. "
503 "The output arrays must each have length at least " << numEntInLclRow
504 <<
" for local row " << LocalRow <<
" on Process "
505 << localRowMap_->getComm ()->getRank () <<
".");
507 else if (numEntInLclRow == static_cast<size_t> (0)) {
526 size_t numEntInMat = 0;
527 A_->getLocalRowCopy (LocalRow, localIndices_ (), Values_ (), numEntInMat);
532 const map_type& matrixDomMap = * (A_->getDomainMap ());
533 const map_type& matrixColMap = * (A_->getColMap ());
535 const size_t capacity =
static_cast<size_t> (std::min (Indices.
size (),
538 const size_t numRows = localRowMap_->getNodeNumElements ();
539 const bool buggy =
true;
540 for (
size_t j = 0; j < numEntInMat; ++j) {
546 const LO matrixLclCol = localIndices_[j];
547 const GO gblCol = matrixColMap.getGlobalElement (matrixLclCol);
556 if ((
size_t) localIndices_[j] < numRows) {
557 Indices[NumEntries] = localIndices_[j];
558 Values[NumEntries] = Values_[j];
562 if (matrixDomMap.isNodeGlobalElement (gblCol)) {
566 if (NumEntries < capacity) {
567 Indices[NumEntries] = matrixLclCol;
568 Values[NumEntries] = Values_[j];
577 template<
class MatrixType>
585 "Ifpack2::LocalFilter does not implement getGlobalRowView.");
589 template<
class MatrixType>
597 "Ifpack2::LocalFilter does not implement getLocalRowView.");
601 template<
class MatrixType>
604 getLocalDiagCopy (Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& diag)
const
608 global_ordinal_type, node_type> vector_type;
611 RCP<vector_type> diagView = diag.offsetViewNonConst (A_->getRowMap (), 0);
612 A_->getLocalDiagCopy (*diagView);
616 template<
class MatrixType>
619 leftScale (
const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& )
622 "Ifpack2::LocalFilter does not implement leftScale.");
626 template<
class MatrixType>
629 rightScale (
const Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& )
632 "Ifpack2::LocalFilter does not implement rightScale.");
636 template<
class MatrixType>
639 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &X,
640 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &Y,
645 typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
647 X.getNumVectors() != Y.getNumVectors(), std::runtime_error,
648 "Ifpack2::LocalFilter::apply: X and Y must have the same number of columns. "
649 "X has " << X.getNumVectors () <<
" columns, but Y has "
650 << Y.getNumVectors () <<
" columns.");
652 #ifdef HAVE_IFPACK2_DEBUG
658 for (
size_t j = 0; j < X.getNumVectors (); ++j) {
659 if (STM::isnaninf (norms[j])) {
664 TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error,
"Ifpack2::LocalFilter::apply: The 1-norm of the input X is NaN or Inf.");
666 #endif // HAVE_IFPACK2_DEBUG
675 applyNonAliased (X_copy, Y, mode, alpha, beta);
677 applyNonAliased (X, Y, mode, alpha, beta);
680 #ifdef HAVE_IFPACK2_DEBUG
686 for (
size_t j = 0; j < Y.getNumVectors (); ++j) {
687 if (STM::isnaninf (norms[j])) {
692 TEUCHOS_TEST_FOR_EXCEPTION( ! good, std::runtime_error,
"Ifpack2::LocalFilter::apply: The 1-norm of the output Y is NaN or Inf.");
694 #endif // HAVE_IFPACK2_DEBUG
697 template<
class MatrixType>
700 applyNonAliased (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &X,
701 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &Y,
704 scalar_type beta)
const
710 const scalar_type zero = STS::zero ();
711 const scalar_type one = STS::one ();
716 else if (beta != one) {
720 const size_t NumVectors = Y.getNumVectors ();
721 const size_t numRows = localRowMap_->getNodeNumElements ();
728 const bool constantStride = X.isConstantStride () && Y.isConstantStride ();
729 if (constantStride) {
732 const size_t x_stride = X.getStride();
733 const size_t y_stride = Y.getStride();
734 ArrayRCP<scalar_type> y_rcp = Y.get1dViewNonConst();
735 ArrayRCP<const scalar_type> x_rcp = X.get1dView();
736 ArrayView<scalar_type> y_ptr = y_rcp();
737 ArrayView<const scalar_type> x_ptr = x_rcp();
738 for (
size_t i = 0; i < numRows; ++i) {
741 getLocalRowCopy (i, localIndices_ (), Values_ (), Nnz);
743 for (
size_t j = 0; j < Nnz; ++j) {
744 const local_ordinal_type col = localIndices_[j];
745 for (
size_t k = 0; k < NumVectors; ++k) {
746 y_ptr[i + y_stride*k] +=
747 alpha * Values_[j] * x_ptr[col + x_stride*k];
752 for (
size_t j = 0; j < Nnz; ++j) {
753 const local_ordinal_type col = localIndices_[j];
754 for (
size_t k = 0; k < NumVectors; ++k) {
755 y_ptr[col + y_stride*k] +=
756 alpha * Values_[j] * x_ptr[i + x_stride*k];
761 for (
size_t j = 0; j < Nnz; ++j) {
762 const local_ordinal_type col = localIndices_[j];
763 for (
size_t k = 0; k < NumVectors; ++k) {
764 y_ptr[col + y_stride*k] +=
765 alpha * STS::conjugate (Values_[j]) * x_ptr[i + x_stride*k];
774 ArrayRCP<ArrayRCP<const scalar_type> > x_ptr = X.get2dView();
775 ArrayRCP<ArrayRCP<scalar_type> > y_ptr = Y.get2dViewNonConst();
777 for (
size_t i = 0; i < numRows; ++i) {
780 getLocalRowCopy (i, localIndices_ (), Values_ (), Nnz);
782 for (
size_t k = 0; k < NumVectors; ++k) {
783 ArrayView<const scalar_type> x_local = (x_ptr())[k]();
784 ArrayView<scalar_type> y_local = (y_ptr())[k]();
785 for (
size_t j = 0; j < Nnz; ++j) {
786 y_local[i] += alpha * Values_[j] * x_local[localIndices_[j]];
791 for (
size_t k = 0; k < NumVectors; ++k) {
792 ArrayView<const scalar_type> x_local = (x_ptr())[k]();
793 ArrayView<scalar_type> y_local = (y_ptr())[k]();
794 for (
size_t j = 0; j < Nnz; ++j) {
795 y_local[localIndices_[j]] += alpha * Values_[j] * x_local[i];
800 for (
size_t k = 0; k < NumVectors; ++k) {
801 ArrayView<const scalar_type> x_local = (x_ptr())[k]();
802 ArrayView<scalar_type> y_local = (y_ptr())[k]();
803 for (
size_t j = 0; j < Nnz; ++j) {
804 y_local[localIndices_[j]] +=
805 alpha * STS::conjugate (Values_[j]) * x_local[i];
813 template<
class MatrixType>
820 template<
class MatrixType>
827 template<
class MatrixType>
829 LocalFilter<MatrixType>::mag_type
832 #ifdef TPETRA_HAVE_KOKKOS_REFACTOR
833 typedef Kokkos::Details::ArithTraits<scalar_type> STS;
834 typedef Kokkos::Details::ArithTraits<mag_type> STM;
841 const size_type maxNumRowEnt = getNodeMaxNumRowEntries ();
844 const size_t numRows =
static_cast<size_t> (localRowMap_->getNodeNumElements ());
847 mag_type sumSquared = STM::zero ();
848 for (
size_t i = 0; i < numRows; ++i) {
849 size_t numEntries = 0;
850 this->getLocalRowCopy (i, ind (), val (), numEntries);
851 for (size_type k = 0; k < static_cast<size_type> (numEntries); ++k) {
852 const mag_type v_k_abs = STS::magnitude (val[k]);
853 sumSquared += v_k_abs * v_k_abs;
856 return STM::squareroot (sumSquared);
859 template<
class MatrixType>
864 std::ostringstream os;
866 os <<
"Ifpack2::LocalFilter: {";
867 os <<
"MatrixType: " << TypeNameTraits<MatrixType>::name ();
868 if (this->getObjectLabel () !=
"") {
869 os <<
", Label: \"" << this->getObjectLabel () <<
"\"";
871 os <<
", Number of rows: " << getGlobalNumRows ()
872 <<
", Number of columns: " << getGlobalNumCols ()
878 template<
class MatrixType>
895 out <<
"Ifpack2::LocalFilter:" << endl;
897 out <<
"MatrixType: " << TypeNameTraits<MatrixType>::name () << endl;
898 if (this->getObjectLabel () !=
"") {
899 out <<
"Label: \"" << this->getObjectLabel () <<
"\"" << endl;
901 out <<
"Number of rows: " << getGlobalNumRows () << endl
902 <<
"Number of columns: " << getGlobalNumCols () << endl
903 <<
"Number of nonzeros: " << NumNonzeros_ << endl;
906 out <<
"Row Map:" << endl;
907 localRowMap_->describe (out, vl);
908 out <<
"Domain Map:" << endl;
909 localDomainMap_->describe (out, vl);
910 out <<
"Range Map:" << endl;
911 localRangeMap_->describe (out, vl);
916 template<
class MatrixType>
917 Teuchos::RCP<
const Tpetra::RowMatrix<
typename MatrixType::scalar_type,
918 typename MatrixType::local_ordinal_type,
919 typename MatrixType::global_ordinal_type,
920 typename MatrixType::node_type> >
929 #define IFPACK2_LOCALFILTER_INSTANT(S,LO,GO,N) \
930 template class Ifpack2::LocalFilter< Tpetra::RowMatrix<S, LO, GO, N> >;
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_LocalFilter_def.hpp:629
virtual Teuchos::RCP< const map_type > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
Definition: Ifpack2_LocalFilter_def.hpp:244
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object to the given output stream.
Definition: Ifpack2_LocalFilter_def.hpp:881
basic_OSTab< char > OSTab
virtual size_t getNodeNumRows() const
The number of rows owned on the calling process.
Definition: Ifpack2_LocalFilter_def.hpp:308
virtual Teuchos::RCP< const map_type > getDomainMap() const
Returns the Map that describes the domain distribution in this matrix.
Definition: Ifpack2_LocalFilter_def.hpp:264
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual bool isFillComplete() const
Returns true if fillComplete() has been called.
Definition: Ifpack2_LocalFilter_def.hpp:421
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Returns the communicator.
Definition: Ifpack2_LocalFilter_def.hpp:232
virtual size_t getNodeMaxNumRowEntries() const
The maximum number of entries across all rows/columns on this process.
Definition: Ifpack2_LocalFilter_def.hpp:393
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
Compute Y = beta*Y + alpha*A_local*X.
Definition: Ifpack2_LocalFilter_def.hpp:639
virtual bool hasTransposeApply() const
Whether this operator supports applying the transpose or conjugate transpose.
Definition: Ifpack2_LocalFilter_def.hpp:814
virtual bool isLocallyIndexed() const
Whether the underlying sparse matrix is locally (opposite of globally) indexed.
Definition: Ifpack2_LocalFilter_def.hpp:407
virtual void getGlobalRowView(global_ordinal_type GlobalRow, Teuchos::ArrayView< const global_ordinal_type > &indices, Teuchos::ArrayView< const scalar_type > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
Definition: Ifpack2_LocalFilter_def.hpp:580
virtual void getLocalDiagCopy(Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &diag) const
Get the diagonal entries of the (locally filtered) matrix.
Definition: Ifpack2_LocalFilter_def.hpp:604
Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > map_type
Type of the Tpetra::Map specialization that this class uses.
Definition: Ifpack2_LocalFilter_decl.hpp:204
virtual mag_type getFrobeniusNorm() const
The Frobenius norm of the (locally filtered) matrix.
Definition: Ifpack2_LocalFilter_def.hpp:830
virtual Teuchos::RCP< const Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > > getGraph() const
The (locally filtered) matrix's graph.
Definition: Ifpack2_LocalFilter_def.hpp:284
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_LocalFilter_decl.hpp:184
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual bool supportsRowViews() const
Returns true if RowViews are supported.
Definition: Ifpack2_LocalFilter_def.hpp:821
virtual size_t getNumEntriesInGlobalRow(global_ordinal_type globalRow) const
The current number of entries on this node in the specified global row.
Definition: Ifpack2_LocalFilter_def.hpp:346
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_LocalFilter_def.hpp:619
virtual global_ordinal_type getIndexBase() const
Returns the index base for global indices for this matrix.
Definition: Ifpack2_LocalFilter_def.hpp:323
virtual bool isGloballyIndexed() const
Whether the underlying sparse matrix is globally (opposite of locally) indexed.
Definition: Ifpack2_LocalFilter_def.hpp:414
virtual size_t getGlobalMaxNumRowEntries() const
The maximum number of entries across all rows/columns on all processes.
Definition: Ifpack2_LocalFilter_def.hpp:386
virtual void getLocalRowCopy(local_ordinal_type LocalRow, const Teuchos::ArrayView< local_ordinal_type > &Indices, const Teuchos::ArrayView< scalar_type > &Values, size_t &NumEntries) const
Get the entries in the given row, using local indices.
Definition: Ifpack2_LocalFilter_def.hpp:474
void resize(size_type new_size, const value_type &x=value_type())
virtual size_t getNodeNumCols() const
The number of columns in the (locally filtered) matrix.
Definition: Ifpack2_LocalFilter_def.hpp:315
LocalFilter(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition: Ifpack2_LocalFilter_def.hpp:88
virtual ~LocalFilter()
Destructor.
Definition: Ifpack2_LocalFilter_def.hpp:226
virtual global_size_t getGlobalNumRows() const
The number of global rows in this matrix.
Definition: Ifpack2_LocalFilter_def.hpp:294
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_LocalFilter_decl.hpp:181
virtual std::string description() const
A one-line description of this object.
Definition: Ifpack2_LocalFilter_def.hpp:861
virtual size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
Definition: Ifpack2_LocalFilter_def.hpp:337
virtual size_t getNumEntriesInLocalRow(local_ordinal_type localRow) const
The current number of entries on this node in the specified local row.
Definition: Ifpack2_LocalFilter_def.hpp:365
virtual Teuchos::RCP< const map_type > getRangeMap() const
Returns the Map that describes the range distribution in this matrix.
Definition: Ifpack2_LocalFilter_def.hpp:274
virtual global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
Definition: Ifpack2_LocalFilter_def.hpp:330
virtual Teuchos::RCP< const row_matrix_type > getUnderlyingMatrix() const
Return matrix that LocalFilter was built on.
Definition: Ifpack2_LocalFilter_def.hpp:921
virtual void getGlobalRowCopy(global_ordinal_type GlobalRow, const Teuchos::ArrayView< global_ordinal_type > &Indices, const Teuchos::ArrayView< scalar_type > &Values, size_t &NumEntries) const
Get the entries in the given row, using global indices.
Definition: Ifpack2_LocalFilter_def.hpp:430
Access only local rows and columns of a sparse matrix.
Definition: Ifpack2_LocalFilter_decl.hpp:160
virtual global_size_t getGlobalNumCols() const
The number of global columns in this matrix.
Definition: Ifpack2_LocalFilter_def.hpp:301
virtual Teuchos::RCP< const map_type > getColMap() const
Returns the Map that describes the column distribution in this matrix.
Definition: Ifpack2_LocalFilter_def.hpp:254
virtual void getLocalRowView(local_ordinal_type LocalRow, Teuchos::ArrayView< const local_ordinal_type > &indices, Teuchos::ArrayView< const scalar_type > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
Definition: Ifpack2_LocalFilter_def.hpp:592
virtual bool hasColMap() const
Whether this matrix has a well-defined column Map.
Definition: Ifpack2_LocalFilter_def.hpp:400