43#ifndef IFPACK2_CONTAINER_DEF_HPP
44#define IFPACK2_CONTAINER_DEF_HPP
47#include <Teuchos_Time.hpp>
53template<
class MatrixType>
55 const Teuchos::RCP<const row_matrix_type>& matrix,
56 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
68 using Teuchos::ArrayView;
85 #ifdef HAVE_IFPACK2_DEBUG
90 Teuchos::ArrayView<const LO> blockRows =
getBlockRows(i);
93 LO row = blockRows[j];
99 TEUCHOS_TEST_FOR_EXCEPTION(
100 !rowMap.isNodeLocalElement(row),
101 std::invalid_argument,
"Ifpack2::Container: "
102 "On process " << rowMap.getComm()->getRank() <<
" of "
103 << rowMap.getComm()->getSize() <<
", in the given set of local row "
104 "indices blockRows = " << Teuchos::toString(blockRows) <<
", the following "
105 "entries is not valid local row index on the calling process: "
112template<
class MatrixType>
116template<
class MatrixType>
117Teuchos::ArrayView<const typename MatrixType::local_ordinal_type>
120 return Teuchos::ArrayView<const LO>
124template<
class MatrixType>
129 LO totalBlockRows = 0;
136 LO rowsInBlock = partitions[i].size();
139 totalBlockRows += rowsInBlock;
140 maxBlockSize_ = std::max(maxBlockSize_, rowsInBlock *
scalarsPerRow_);
154template<
class MatrixType>
155void Container<MatrixType>::getMatDiag()
const
164template<
class MatrixType>
169template<
class MatrixType>
174template<
class MatrixType>
176applyMV(
const mv_type& X, mv_type& Y)
const
178 TEUCHOS_TEST_FOR_EXCEPT_MSG(
true,
"Not implemented.");
181template<
class MatrixType>
185 TEUCHOS_TEST_FOR_EXCEPT_MSG(
true,
"Not implemented.");
188template<
class MatrixType>
195template<
class MatrixType>
197 SC dampingFactor, LO i)
const
199 TEUCHOS_TEST_FOR_EXCEPT_MSG(
true,
"Not implemented.");
202template <
class MatrixType>
205 using STS = Teuchos::ScalarTraits<ISC>;
206 const ISC one = STS::one();
208 size_t numVecs = X.extent(1);
210 for (LO i = 0; i < numBlocks_; i++)
213 if(blockSizes_[i] != 1 || hasBlockCrs_)
215 if(blockSizes_[i] == 0 )
217 apply(X, Y, i, Teuchos::NO_TRANS, dampingFactor, one);
221 LO LRID = blockRows_[blockOffsets_[i]];
223 auto diagView = Diag_->getLocalViewHost(Tpetra::Access::ReadOnly);
224 ISC d = one * (
static_cast<ISC
> (dampingFactor)) / diagView(LRID, 0);
225 for(
size_t nv = 0; nv < numVecs; nv++)
228 Y(LRID, nv) += x * d;
234template <
class MatrixType>
237 using STS = Teuchos::ScalarTraits<SC>;
239 for(LO i = 0; i < numBlocks_; i++)
242 if(blockSizes_[i] == 0)
244 if(blockSizes_[i] != 1) {
246 weightedApply(X, Y, W, i, Teuchos::NO_TRANS, dampingFactor, STS::one());
253 HostView tempo(
"", X.extent(0), X.extent(1));
254 size_t numVecs = X.extent(1);
255 LO bOffset = blockOffsets_[i];
256 for (LO ii = 0; ii < blockSizes_[i]; ii++) {
257 LO LRID = blockRows_[bOffset++];
258 for (
size_t jj = 0; jj < numVecs; jj++) tempo(LRID,jj)=X(LRID,jj)/ W(LRID,0);
260 weightedApply(tempo, Y, W, i, Teuchos::NO_TRANS, dampingFactor, STS::one());
265 const ISC one = STS::one();
266 size_t numVecs = X.extent(1);
269 auto diagView =
Diag_->getLocalViewHost(Tpetra::Access::ReadOnly);
270 ISC recip = one / diagView(LRID, 0);
271 ISC wval = W(LRID,0);
272 ISC combo = wval*recip;
273 ISC d = combo*(
static_cast<ISC> (dampingFactor));
274 for(
size_t nv = 0; nv < numVecs; nv++)
277 Y(LRID, nv) = x * d + Y(LRID, nv);
285template<
class MatrixType,
typename LocalScalarType>
288 SC dampingFactor, LO i)
const
290 using Teuchos::ArrayView;
291 using STS = Teuchos::ScalarTraits<ISC>;
292 size_t numVecs = X.extent(1);
293 const ISC one = STS::one();
294 if(this->blockSizes_[i] == 0)
300 const size_t localNumRows = this->blockSizes_[i];
301 using inds_type =
typename block_crs_matrix_type::local_inds_host_view_type;
302 using vals_type =
typename block_crs_matrix_type::values_host_view_type;
303 for(
size_t j = 0; j < localNumRows; j++)
305 LO row = blockRows[j];
309 LO numEntries = (LO) colinds.size();
310 for(
size_t m = 0; m < numVecs; m++)
313 Resid(row * this->bcrsBlockSize_ + localR, m) = X(row * this->bcrsBlockSize_ + localR, m);
314 for (LO k = 0; k < numEntries; ++k)
316 const LO col = colinds[k];
317 for(
int localR = 0; localR < this->bcrsBlockSize_; localR++)
319 for(
int localC = 0; localC < this->bcrsBlockSize_; localC++)
321 Resid(row * this->bcrsBlockSize_ + localR, m) -=
322 values[k * this->bcrsBlockSize_ * this->bcrsBlockSize_ + localR + localC * this->bcrsBlockSize_]
323 * Y2(col * this->bcrsBlockSize_ + localC, m); }
334 apply(Resid, Y2, i, Teuchos::NO_TRANS, dampingFactor, one);
336 else if(!this->
hasBlockCrs_ && this->blockSizes_[i] == 1)
343 using container_exec_space =
typename ContainerImpl<MatrixType, LocalScalarType>::crs_matrix_type::execution_space;
344 container_exec_space().fence();
346 using size_type =
typename crs_matrix_type::local_matrix_host_type::size_type;
347 const auto& rowmap = localA.graph.row_map;
348 const auto& entries = localA.graph.entries;
349 const auto& values = localA.values;
351 auto diagView = this->
Diag_->getLocalViewHost(Tpetra::Access::ReadOnly);
352 ISC d = (
static_cast<ISC> (dampingFactor)) / diagView(LRID, 0);
353 for(
size_t m = 0; m < numVecs; m++)
357 for(size_type k = rowmap(LRID); k < rowmap(LRID + 1); k++) {
358 const LO col = entries(k);
359 r -= values(k) * Y2(col, m);
367 std::is_same<typename crs_matrix_type::device_type::memory_space, Kokkos::HostSpace>::value)
371 using container_exec_space =
typename ContainerImpl<MatrixType, LocalScalarType>::crs_matrix_type::execution_space;
372 container_exec_space().fence();
374 using size_type =
typename crs_matrix_type::local_matrix_host_type::size_type;
375 const auto& rowmap = localA.graph.row_map;
376 const auto& entries = localA.graph.entries;
377 const auto& values = localA.values;
379 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
381 const LO row = blockRows[j];
382 for(
size_t m = 0; m < numVecs; m++)
385 for(size_type k = rowmap(row); k < rowmap(row + 1); k++)
387 const LO col = entries(k);
388 r -= values(k) * Y2(col, m);
399 apply(Resid, Y2, i, Teuchos::NO_TRANS, dampingFactor, one);
406 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
408 const LO row = blockRows[j];
410 for(
size_t m = 0; m < numVecs; m++)
412 Resid(row, m) = X(row, m);
413 for (
size_t k = 0; k < rowView.size(); ++k)
415 const LO col = rowView.ind(k);
416 Resid(row, m) -= rowView.val(k) * Y2(col, m);
426 apply(Resid, Y2, i, Teuchos::NO_TRANS, dampingFactor, one);
430template<
class MatrixType>
431void Container<MatrixType>::
432DoGaussSeidel(ConstHostView X, HostView Y, HostView Y2, SC dampingFactor)
const
434 using Teuchos::Array;
435 using Teuchos::ArrayRCP;
436 using Teuchos::ArrayView;
440 using Teuchos::rcpFromRef;
443 auto numVecs = X.extent(1);
445 HostView Resid(
"", X.extent(0), X.extent(1));
448 DoGSBlock(X, Y, Y2, Resid, dampingFactor, i);
453 for (
size_t m = 0; m < numVecs; ++m)
463template<
class MatrixType>
464void Container<MatrixType>::
465DoSGS(ConstHostView X, HostView Y, HostView Y2, SC dampingFactor)
const
468 using Teuchos::Array;
469 using Teuchos::ArrayRCP;
470 using Teuchos::ArrayView;
475 using Teuchos::rcpFromRef;
476 auto numVecs = X.extent(1);
477 HostView Resid(
"", X.extent(0), X.extent(1));
481 DoGSBlock(X, Y, Y2, Resid, dampingFactor, i);
483 static_assert(std::is_signed<LO>::value,
484 "Local ordinal must be signed (unsigned breaks reverse iteration to 0)");
488 DoGSBlock(X, Y, Y2, Resid, dampingFactor, i);
493 for (
size_t m = 0; m < numVecs; ++m)
503template<
class MatrixType>
510 blockOffsets_.clear();
511 Diag_ = Teuchos::null;
516template<
class MatrixType,
class LocalScalarType>
519 const Teuchos::RCP<const row_matrix_type>& matrix,
520 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
522 :
Container<MatrixType>(matrix, partitions, pointIndexed) {}
524template<
class MatrixType,
class LocalScalarType>
528template<
class MatrixType,
class LocalScalarType>
532template<
class MatrixType,
class LocalScalarType>
539 TEUCHOS_TEST_FOR_EXCEPT_MSG(
true,
"Not implemented.");
542template<
class MatrixType,
class LocalScalarType>
544applyMV (
const mv_type& X, mv_type& Y)
const
547 HostView YView = Y.getLocalViewHost(Tpetra::Access::ReadWrite);
548 this->apply (XView, YView, 0);
551template<
class MatrixType,
class LocalScalarType>
557 ConstHostView XView = X.getLocalViewHost(Tpetra::Access::ReadOnly);
558 HostView YView = Y.getLocalViewHost(Tpetra::Access::ReadWrite);
559 ConstHostView WView = W.getLocalViewHost(Tpetra::Access::ReadOnly);
560 weightedApply (XView, YView, WView, 0);
563template<
class MatrixType,
class LocalScalarType>
570template<
class MatrixType,
class LocalScalarType>
575 Teuchos::ETransp mode,
577 const LSC beta)
const
579 TEUCHOS_TEST_FOR_EXCEPT_MSG(
true,
"Not implemented.");
582template<
class MatrixType,
class LocalScalarType>
583typename ContainerImpl<MatrixType, LocalScalarType>::LO
587 LO LO_INVALID = Teuchos::OrdinalTraits<LO>::invalid();
588 GO GO_INVALID = Teuchos::OrdinalTraits<GO>::invalid();
589 const map_type& globalRowMap = *(this->
inputMatrix_->getRowMap());
590 const map_type& globalColMap = *(this->
inputMatrix_->getColMap());
598 GO diagGID = globalRowMap.getGlobalElement(rowLID);
599 TEUCHOS_TEST_FOR_EXCEPTION(
600 diagGID == GO_INVALID,
601 std::runtime_error,
"Ifpack2::Container::translateRowToCol: "
602 "On process " << this->
inputMatrix_->getRowMap()->getComm()->getRank() <<
603 ", at least one row index in the set of local "
604 "row indices given to the constructor is not a valid local row index in "
605 "the input matrix's row Map on this process. This should be impossible "
606 "because the constructor checks for this case. Here is the complete set "
607 "of invalid local row indices: " << rowLID <<
". "
608 "Please report this bug to the Ifpack2 developers.");
610 LO colLID = globalColMap.getLocalElement(diagGID);
611 TEUCHOS_TEST_FOR_EXCEPTION(
612 colLID == LO_INVALID,
613 std::runtime_error,
"Ifpack2::Container::translateRowToCol: "
614 "On process " << this->
inputMatrix_->getRowMap()->getComm()->getRank() <<
", "
615 "at least one row index in the set of row indices given to the constructor "
616 "does not have a corresponding column index in the input matrix's column "
617 "Map. This probably means that the column(s) in question is/are empty on "
618 "this process, which would make the submatrix to extract structurally "
619 "singular. The invalid global column index is " << diagGID <<
".");
626template<
class MatrixType,
class LocalScalarType>
628apply (ConstHostView X,
631 Teuchos::ETransp mode,
635 using Teuchos::ArrayView;
648 TEUCHOS_TEST_FOR_EXCEPTION(
649 ! this->
IsComputed_, std::runtime_error,
"Ifpack2::Container::apply: "
650 "You must have called the compute() method before you may call apply(). "
651 "You may call the apply() method as many times as you want after calling "
652 "compute() once, but you must have called compute() at least once.");
654 const size_t numVecs = X.extent(1);
693 X_local_ = HostViewLocal(
"X_local", this->
blockRows_.size() * this->scalarsPerRow_, numVecs);
694 Y_local_ = HostViewLocal(
"Y_local", this->
blockRows_.size() * this->scalarsPerRow_, numVecs);
698 for(
int i = 0; i < this->numBlocks_; i++)
701 (this->
blockOffsets_[i] + this->blockSizes_[i]) * this->scalarsPerRow_);
703 Y_localBlocks_.emplace_back(Y_local_, blockBounds, Kokkos::ALL());
707 const ArrayView<const LO> blockRows = this->
getBlockRows(blockIndex);
732 mvgs.scatterViewToView (Y,
Y_localBlocks_[blockIndex], blockRows);
737template<
class MatrixType,
class LocalScalarType>
743 Teuchos::ETransp mode,
747 using Teuchos::ArrayRCP;
748 using Teuchos::ArrayView;
749 using Teuchos::Range1D;
754 using Teuchos::rcp_const_cast;
756 using STS = Teuchos::ScalarTraits<SC>;
765 const char prefix[] =
"Ifpack2::Container::weightedApply: ";
766 TEUCHOS_TEST_FOR_EXCEPTION(
767 ! this->
IsComputed_, std::runtime_error, prefix <<
"You must have called the "
768 "compute() method before you may call this method. You may call "
769 "weightedApply() as many times as you want after calling compute() once, "
770 "but you must have called compute() at least once first.");
773 TEUCHOS_TEST_FOR_EXCEPTION(
774 this->
scalarsPerRow_ > 1, std::logic_error, prefix <<
"Use of block rows isn't allowed "
775 "in overlapping Jacobi (the only method that uses weightedApply");
777 const size_t numVecs = X.extent(1);
779 TEUCHOS_TEST_FOR_EXCEPTION(
780 X.extent(1) != Y.extent(1), std::runtime_error,
781 prefix <<
"X and Y have different numbers of vectors (columns). X has "
782 << X.extent(1) <<
", but Y has " << Y.extent(1) <<
".");
788 const size_t numRows = this->blockSizes_[blockIndex];
821 X_local_ = HostViewLocal(
"X_local", this->
blockRows_.size() * this->scalarsPerRow_, numVecs);
822 Y_local_ = HostViewLocal(
"Y_local", this->
blockRows_.size() * this->scalarsPerRow_, numVecs);
826 for(
int i = 0; i < this->numBlocks_; i++)
829 (this->
blockOffsets_[i] + this->blockSizes_[i]) * this->scalarsPerRow_);
831 Y_localBlocks_.emplace_back(Y_local_, blockBounds, Kokkos::ALL());
840 ArrayView<const LO> blockRows = this->
getBlockRows(blockIndex);
862 auto maxBS = this->maxBlockSize_;
866 mvgs.gatherViewToView (D_local, D, blockRows);
868 for(
size_t j = 0; j < numVecs; j++)
869 for(
size_t i = 0; i < numRows; i++)
870 X_scaled(i, j) =
X_localBlocks_[blockIndex](i, j) * D_local(i, 0);
872 HostSubviewLocal Y_temp(
weightedApplyScratch_, std::make_pair(maxBS * 2, maxBS * 2 + bs), Kokkos::ALL());
874 this->
solveBlock (X_scaled, Y_temp, blockIndex, mode, STS::one(), STS::zero());
882 for(
size_t j = 0; j < numVecs; j++)
883 for(
size_t i = 0; i < numRows; i++)
888 mvgs.scatterViewToView (Y,
Y_localBlocks_[blockIndex], blockRows);
891template<
class MatrixType,
class LocalScalarType>
893 typename ContainerImpl<MatrixType, LocalScalarType>::SC,
894 typename ContainerImpl<MatrixType, LocalScalarType>::LO,
895 typename ContainerImpl<MatrixType, LocalScalarType>::GO,
896 typename ContainerImpl<MatrixType, LocalScalarType>::NO>
901 typedef typename MatrixType::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
902 typedef typename MatrixType::nonconst_values_host_view_type nonconst_values_host_view_type;
904 typedef typename MatrixType::local_inds_host_view_type local_inds_host_view_type;
905 typedef typename MatrixType::values_host_view_type values_host_view_type;
906 using IST =
typename row_matrix_type::impl_scalar_type;
910 using h_inds_type =
typename block_crs_matrix_type::local_inds_host_view_type;
911 using h_vals_type =
typename block_crs_matrix_type::values_host_view_type;
915 size_t numEntries = colinds.size();
918 h_vals_type subvals = Kokkos::subview(values,std::pair<size_t,size_t>(row % this->
bcrsBlockSize_,values.size()));
919 return StridedRowView(subvals, colinds, this->bcrsBlockSize_, numEntries * this->bcrsBlockSize_);
923 size_t maxEntries = this->
inputMatrix_->getLocalMaxNumRowEntries();
924 Teuchos::Array<LO> inds(maxEntries);
925 Teuchos::Array<SC> vals(maxEntries);
926 nonconst_local_inds_host_view_type inds_v(inds.data(),maxEntries);
927 nonconst_values_host_view_type vals_v(
reinterpret_cast<IST*
>(vals.data()),maxEntries);
929 this->
inputMatrix_->getLocalRowCopy(row, inds_v, vals_v, numEntries);
930 vals.resize(numEntries); inds.resize(numEntries);
931 return StridedRowView(vals, inds);
936 local_inds_host_view_type colinds;
937 values_host_view_type values;
938 this->
inputMatrix_->getLocalRowView(row, colinds, values);
939 return StridedRowView(values, colinds, 1, colinds.size());
943template<
class MatrixType,
class LocalScalarType>
947 X_localBlocks_.clear();
948 Y_localBlocks_.clear();
949 X_local_ = HostViewLocal();
950 Y_local_ = HostViewLocal();
951 Container<MatrixType>::clearBlocks();
957template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
959StridedRowView(h_vals_type vals_, h_inds_type inds_,
int blockSize_,
size_t nnz_)
960 : vals(vals_), inds(inds_), blockSize(blockSize_), nnz(nnz_)
963template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
965StridedRowView(Teuchos::Array<SC>& vals_, Teuchos::Array<LO>& inds_)
966 : vals(), inds(), blockSize(1), nnz(vals_.size())
968 valsCopy.swap(vals_);
969 indsCopy.swap(inds_);
972template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
976 #ifdef HAVE_IFPACK2_DEBUG
977 TEUCHOS_TEST_FOR_EXCEPTION(i >= nnz, std::runtime_error,
978 "Out-of-bounds access into Ifpack2::Container::StridedRowView");
985 return vals[i * blockSize];
991template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
992LocalOrdinal StridedRowView<Scalar, LocalOrdinal, GlobalOrdinal, Node>::
995 #ifdef HAVE_IFPACK2_DEBUG
996 TEUCHOS_TEST_FOR_EXCEPTION(i >= nnz, std::runtime_error,
997 "Out-of-bounds access into Ifpack2::Container::StridedRowView");
1005 return inds[i / blockSize] * blockSize + i % blockSize;
1011template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1021template <
class MatrixType>
1024 return obj.
print(os);
1027#define IFPACK2_CONTAINER_INSTANT(S,LO,GO,N) \
1028 template class Ifpack2::Container<Tpetra::RowMatrix<S, LO, GO, N>>; \
1029 template class Ifpack2::ContainerImpl<Tpetra::RowMatrix<S, LO, GO, N>, S>; \
1030 template class Ifpack2::Details::StridedRowView<S, LO, GO, N>; \
1031 template std::ostream& operator<< <Tpetra::RowMatrix<S, LO, GO, N>>( \
1032 std::ostream& os, const Ifpack2::Container<Tpetra::RowMatrix<S, LO, GO, N>>& obj);
Declaration and definition of the Ifpack2::Details::MultiVectorLocalGatherScatter class.
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS)....
Definition Ifpack2_Container_decl.hpp:344
std::vector< HostSubviewLocal > X_localBlocks_
Views for holding pieces of X corresponding to each block.
Definition Ifpack2_Container_decl.hpp:523
HostViewLocal X_local_
Scratch vectors used in apply().
Definition Ifpack2_Container_decl.hpp:507
virtual void setParameters(const Teuchos::ParameterList &List)
Set parameters, if any.
Definition Ifpack2_Container_def.hpp:530
Details::StridedRowView< SC, LO, GO, NO > getInputRowView(LO row) const
Definition Ifpack2_Container_def.hpp:898
std::vector< HostSubviewLocal > Y_localBlocks_
Views for holding pieces of Y corresponding to each block.
Definition Ifpack2_Container_decl.hpp:526
virtual void weightedApply(ConstHostView X, HostView Y, ConstHostView D, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, SC alpha=Teuchos::ScalarTraits< SC >::one(), SC beta=Teuchos::ScalarTraits< SC >::zero()) const
Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y.
Definition Ifpack2_Container_def.hpp:739
void applyMV(const mv_type &X, mv_type &Y) const
Wrapper for apply with MVs, used in unit tests (never called by BlockRelaxation)
Definition Ifpack2_Container_def.hpp:544
virtual void apply(ConstHostView X, HostView Y, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, SC alpha=Teuchos::ScalarTraits< SC >::one(), SC beta=Teuchos::ScalarTraits< SC >::zero()) const
Compute Y := alpha * M^{-1} X + beta*Y.
Definition Ifpack2_Container_def.hpp:628
LocalScalarType LSC
The internal representation of LocalScalarType in Kokkos::View.
Definition Ifpack2_Container_decl.hpp:363
void DoGSBlock(ConstHostView X, HostView Y, HostView Y2, HostView Resid, SC dampingFactor, LO i) const
Do one step of Gauss-Seidel on block i (used by DoGaussSeidel and DoSGS)
Definition Ifpack2_Container_def.hpp:286
virtual void solveBlock(ConstHostSubviewLocal X, HostSubviewLocal Y, int blockIndex, Teuchos::ETransp mode, const LSC alpha, const LSC beta) const
Definition Ifpack2_Container_def.hpp:572
HostViewLocal weightedApplyScratch_
Definition Ifpack2_Container_decl.hpp:520
static std::string getName()
Definition Ifpack2_Container_def.hpp:565
LO translateRowToCol(LO row)
Definition Ifpack2_Container_def.hpp:585
virtual ~ContainerImpl()
Destructor.
Definition Ifpack2_Container_def.hpp:526
void weightedApplyMV(const mv_type &X, mv_type &Y, vector_type &W) const
Wrapper for weightedApply with MVs, used in unit tests (never called by BlockRelaxation)
Definition Ifpack2_Container_def.hpp:553
virtual void applyInverseJacobi(const mv_type &, mv_type &, SC dampingFactor, bool, int) const
Compute Y := (1 - a) Y + a D^{-1} (X - R*Y).
Definition Ifpack2_Container_def.hpp:534
virtual ~Container()
Destructor.
Definition Ifpack2_Container_def.hpp:114
bool hasBlockCrs_
Whether the input matrix is a BlockCRS matrix.
Definition Ifpack2_Container_decl.hpp:310
Teuchos::RCP< const crs_matrix_type > inputCrsMatrix_
The input matrix, dynamic cast to CrsMatrix. May be null.
Definition Ifpack2_Container_decl.hpp:286
bool IsParallel_
Whether the problem is distributed across multiple MPI processes.
Definition Ifpack2_Container_decl.hpp:302
int numBlocks_
The number of blocks (partitions) in the container.
Definition Ifpack2_Container_decl.hpp:292
Teuchos::ArrayView< const LO > getBlockRows(int blockIndex) const
Local indices of the rows of the input matrix that belong to this block.
Definition Ifpack2_Container_def.hpp:118
static std::string getName()
Definition Ifpack2_Container_def.hpp:190
Teuchos::RCP< vector_type > Diag_
Definition Ifpack2_Container_decl.hpp:300
virtual void weightedApply(ConstHostView X, HostView Y, ConstHostView D, int blockIndex, Teuchos::ETransp mode=Teuchos::NO_TRANS, SC alpha=Teuchos::ScalarTraits< SC >::one(), SC beta=Teuchos::ScalarTraits< SC >::zero()) const=0
int bcrsBlockSize_
If hasBlockCrs_, the number of DOFs per vertex. Otherwise 1.
Definition Ifpack2_Container_decl.hpp:312
Teuchos::RCP< const block_crs_matrix_type > inputBlockMatrix_
The input matrix, dynamic cast to BlockCrsMatrix. May be null.
Definition Ifpack2_Container_decl.hpp:289
typename Kokkos::Details::ArithTraits< SC >::val_type ISC
Definition Ifpack2_Container_decl.hpp:135
virtual void DoGSBlock(ConstHostView X, HostView Y, HostView Y2, HostView Resid, SC dampingFactor, LO i) const
Do one step of Gauss-Seidel on block i (used by DoGaussSeidel and DoSGS)
Definition Ifpack2_Container_def.hpp:196
virtual void applyMV(const mv_type &X, mv_type &Y) const
Wrapper for apply with MultiVector.
Definition Ifpack2_Container_def.hpp:176
LO NumLocalRows_
Number of local rows in input matrix.
Definition Ifpack2_Container_decl.hpp:304
bool isInitialized() const
Whether the container has been successfully initialized.
Definition Ifpack2_Container_def.hpp:165
void setBlockSizes(const Teuchos::Array< Teuchos::Array< LO > > &partitions)
Initialize arrays with information about block sizes.
Definition Ifpack2_Container_def.hpp:125
LO scalarsPerRow_
Definition Ifpack2_Container_decl.hpp:317
Teuchos::Array< LO > blockSizes_
Number of rows in each block.
Definition Ifpack2_Container_decl.hpp:296
GO NumGlobalNonzeros_
Number of nonzeros in input matrix.
Definition Ifpack2_Container_decl.hpp:308
Container(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, bool pointIndexed)
Constructor.
Definition Ifpack2_Container_def.hpp:54
virtual std::ostream & print(std::ostream &os) const =0
Print basic information about the container to os.
virtual void weightedApplyMV(const mv_type &X, mv_type &Y, vector_type &W) const
Wrapper for weightedApply with MultiVector.
Definition Ifpack2_Container_def.hpp:183
Teuchos::RCP< const row_matrix_type > inputMatrix_
The input matrix to the constructor.
Definition Ifpack2_Container_decl.hpp:283
bool pointIndexed_
(If hasBlockCrs_) Whether the blocks are described using sub-block row indices instead of full block ...
Definition Ifpack2_Container_decl.hpp:314
bool IsInitialized_
If true, the container has been successfully initialized.
Definition Ifpack2_Container_decl.hpp:320
Teuchos::Array< LO > blockRows_
Local indices of the rows of the input matrix that belong to this block.
Definition Ifpack2_Container_decl.hpp:294
GO NumGlobalRows_
Number of global rows in input matrix.
Definition Ifpack2_Container_decl.hpp:306
bool isComputed() const
Whether the container has been successfully computed.
Definition Ifpack2_Container_def.hpp:170
typename mv_type::dual_view_type::t_host HostView
Definition Ifpack2_Container_decl.hpp:139
Teuchos::Array< LO > blockOffsets_
Starting index in blockRows_ of local row indices for each block.
Definition Ifpack2_Container_decl.hpp:298
bool IsComputed_
If true, the container has been successfully computed.
Definition Ifpack2_Container_decl.hpp:322
Implementation detail of Ifpack2::Container subclasses.
Definition Ifpack2_Details_MultiVectorLocalGatherScatter.hpp:85
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:74
Structure for read-only views of general matrix rows.
Definition Ifpack2_Container_decl.hpp:537
StridedRowView(h_vals_type vals_, h_inds_type inds_, int blockSize_, size_t nnz_)
Constructor for row views (preferred)
Definition Ifpack2_Container_def.hpp:959