10 #ifndef IFPACK2_CONTAINER_DECL_HPP
11 #define IFPACK2_CONTAINER_DECL_HPP
16 #include "Ifpack2_ConfigDefs.hpp"
17 #include "Tpetra_RowMatrix.hpp"
18 #include "Teuchos_Describable.hpp"
19 #include <Tpetra_Map.hpp>
20 #include <Tpetra_BlockCrsMatrix.hpp>
24 #ifndef DOXYGEN_SHOULD_SKIP_THIS
29 #endif // DOXYGEN_SHOULD_SKIP_THIS
78 template<
class MatrixType>
82 using scalar_type =
typename MatrixType::scalar_type;
83 using local_ordinal_type =
typename MatrixType::local_ordinal_type;
84 using global_ordinal_type =
typename MatrixType::global_ordinal_type;
85 using node_type =
typename MatrixType::node_type;
86 using SC = scalar_type;
87 using LO = local_ordinal_type;
88 using GO = global_ordinal_type;
90 using import_type = Tpetra::Import<LO, GO, NO>;
91 using mv_type = Tpetra::MultiVector<SC, LO, GO, NO>;
92 using vector_type = Tpetra::Vector<SC, LO, GO, NO>;
93 using map_type = Tpetra::Map<LO, GO, NO>;
94 using crs_matrix_type = Tpetra::CrsMatrix<SC, LO, GO, NO>;
95 using block_crs_matrix_type = Tpetra::BlockCrsMatrix<SC, LO, GO, NO>;
96 using row_matrix_type = Tpetra::RowMatrix<SC, LO, GO, NO>;
98 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
99 "Ifpack2::Container: Please use MatrixType = Tpetra::RowMatrix.");
102 using ISC =
typename Kokkos::ArithTraits<SC>::val_type;
106 using HostView =
typename mv_type::dual_view_type::t_host;
107 using ConstHostView =
typename HostView::const_type;
157 void getMatDiag()
const;
175 void DoJacobi(ConstHostView X,
HostView Y, SC dampingFactor)
const;
176 void DoOverlappingJacobi(ConstHostView X,
HostView Y, ConstHostView W, SC dampingFactor,
bool nonsymScaling)
const;
177 void DoGaussSeidel(ConstHostView X,
HostView Y,
HostView Y2, SC dampingFactor)
const;
178 void DoSGS(ConstHostView X,
HostView Y,
HostView Y2, SC dampingFactor)
const;
196 apply(ConstHostView X,
227 virtual void applyMV (
const mv_type& X, mv_type& Y)
const;
232 vector_type& W)
const;
234 virtual void clearBlocks();
237 virtual std::ostream&
print (std::ostream& os)
const = 0;
247 SC dampingFactor, LO i)
const;
296 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
309 template<
class MatrixType,
class LocalScalarType>
315 using local_scalar_type = LocalScalarType;
316 using SC =
typename Container<MatrixType>::scalar_type;
317 using LO =
typename Container<MatrixType>::local_ordinal_type;
318 using GO =
typename Container<MatrixType>::global_ordinal_type;
319 using NO =
typename Container<MatrixType>::node_type;
321 using typename Container<MatrixType>::import_type;
322 using typename Container<MatrixType>::row_matrix_type;
323 using typename Container<MatrixType>::crs_matrix_type;
324 using typename Container<MatrixType>::block_crs_matrix_type;
325 using typename Container<MatrixType>::mv_type;
326 using typename Container<MatrixType>::vector_type;
327 using typename Container<MatrixType>::map_type;
330 using LSC = LocalScalarType;
331 using LISC =
typename Kokkos::ArithTraits<LSC>::val_type;
333 using local_mv_type = Tpetra::MultiVector<LSC, LO, GO, NO>;
336 using typename Container<MatrixType>::ConstHostView;
337 using HostViewLocal =
typename local_mv_type::dual_view_type::t_host;
338 using HostSubviewLocal = Kokkos::View<LISC**, Kokkos::LayoutStride, typename HostViewLocal::memory_space>;
339 using ConstHostSubviewLocal = Kokkos::View<const LISC**, Kokkos::LayoutStride, typename HostViewLocal::memory_space>;
341 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
342 "Ifpack2::Container: Please use MatrixType = Tpetra::RowMatrix.");
408 apply(ConstHostView X,
439 void applyMV (
const mv_type& X, mv_type& Y)
const;
444 vector_type& W)
const;
446 virtual void clearBlocks();
449 virtual std::ostream&
print (std::ostream& os)
const = 0;
458 SC dampingFactor, LO i)
const;
471 const LSC beta)
const;
475 mutable HostViewLocal Y_local_;
502 template<
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
503 struct StridedRowView
506 using LO = LocalOrdinal;
508 using block_crs_matrix_type = Tpetra::BlockCrsMatrix<SC, LO, GlobalOrdinal, Node>;
510 using h_inds_type =
typename block_crs_matrix_type::local_inds_host_view_type;
511 using h_vals_type =
typename block_crs_matrix_type::values_host_view_type;
513 StridedRowView(h_vals_type vals_, h_inds_type inds_,
int blockSize_,
size_t nnz_);
521 SC val(
size_t i)
const;
522 LO ind(
size_t i)
const;
540 template <
class MatrixType>
541 std::ostream& operator<<(std::ostream& os, const Ifpack2::Container<MatrixType>& obj);
549 template<
class MatrixType>
553 static std::string name () {
554 return std::string (
"Ifpack2::Container<") +
558 static std::string concreteName (const ::Ifpack2::Container<MatrixType>&) {
565 #endif // IFPACK2_CONTAINER_HPP
virtual void applyInverseJacobi(const mv_type &, mv_type &, SC dampingFactor, bool, int) const =0
Compute Y := (1 - a) Y + a D^{-1} (X - R*Y).
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:595
virtual void initialize()=0
Do all set-up operations that only require matrix structure.
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:520
bool isComputed() const
Whether the container has been successfully computed.
Definition: Ifpack2_Container_def.hpp:137
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:253
std::vector< HostSubviewLocal > X_localBlocks_
Views for holding pieces of X corresponding to each block.
Definition: Ifpack2_Container_decl.hpp:490
Teuchos::Array< LO > blockRows_
Local indices of the rows of the input matrix that belong to this block.
Definition: Ifpack2_Container_decl.hpp:261
void setBlockSizes(const Teuchos::Array< Teuchos::Array< LO > > &partitions)
Initialize arrays with information about block sizes.
Definition: Ifpack2_Container_def.hpp:92
bool IsInitialized_
If true, the container has been successfully initialized.
Definition: Ifpack2_Container_decl.hpp:287
Teuchos::Array< LO > blockSizes_
Number of rows in each block.
Definition: Ifpack2_Container_decl.hpp:263
typename Kokkos::ArithTraits< SC >::val_type ISC
Internal representation of Scalar in Kokkos::View.
Definition: Ifpack2_Container_decl.hpp:102
StridedRowView(h_vals_type vals_, h_inds_type inds_, int blockSize_, size_t nnz_)
Constructor for row views (preferred)
Definition: Ifpack2_Container_def.hpp:926
std::vector< HostSubviewLocal > Y_localBlocks_
Views for holding pieces of Y corresponding to each block.
Definition: Ifpack2_Container_decl.hpp:493
virtual void applyMV(const mv_type &X, mv_type &Y) const
Wrapper for apply with MultiVector.
Definition: Ifpack2_Container_def.hpp:143
virtual void initialize()=0
Do all set-up operations that only require matrix structure.
bool pointIndexed_
(If hasBlockCrs_) Whether the blocks are described using sub-block row indices instead of full block ...
Definition: Ifpack2_Container_decl.hpp:281
int numBlocks_
The number of blocks (partitions) in the container.
Definition: Ifpack2_Container_decl.hpp:259
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS). This class allows a custom scalar type (LocalScalarType) to be used for storing blocks and solving the block systems. Hiding this template parameter from the Container interface simplifies the BlockRelaxation and ContainerFactory classes.
Definition: Ifpack2_Container_decl.hpp:310
virtual void setParameters(const Teuchos::ParameterList &List)
Set parameters, if any.
Definition: Ifpack2_Container_def.hpp:497
bool isInitialized() const
Whether the container has been successfully initialized.
Definition: Ifpack2_Container_def.hpp:132
LO scalarsPerRow_
Definition: Ifpack2_Container_decl.hpp:284
virtual ~Container()
Destructor.
Definition: Ifpack2_Container_def.hpp:81
bool hasBlockCrs_
Whether the input matrix is a BlockCRS matrix.
Definition: Ifpack2_Container_decl.hpp:277
Details::StridedRowView< SC, LO, GO, NO > getInputRowView(LO row) const
View a row of the input matrix.
Definition: Ifpack2_Container_def.hpp:865
HostViewLocal weightedApplyScratch_
Definition: Ifpack2_Container_decl.hpp:487
virtual void compute()=0
Extract the local diagonal blocks and prepare the solver.
InverseType::scalar_type LSC
The internal representation of LocalScalarType in Kokkos::View.
Definition: Ifpack2_Container_decl.hpp:330
Container(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, bool pointIndexed)
Constructor.
Definition: Ifpack2_Container_def.hpp:21
int bcrsBlockSize_
If hasBlockCrs_, the number of DOFs per vertex. Otherwise 1.
Definition: Ifpack2_Container_decl.hpp:279
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 =0
Compute Y := alpha * M^{-1} X + beta*Y.
LO translateRowToCol(LO row)
Definition: Ifpack2_Container_def.hpp:552
HostViewLocal X_local_
Scratch vectors used in apply().
Definition: Ifpack2_Container_decl.hpp:474
Teuchos::RCP< const block_crs_matrix_type > inputBlockMatrix_
The input matrix, dynamic cast to BlockCrsMatrix. May be null.
Definition: Ifpack2_Container_decl.hpp:256
Teuchos::RCP< const row_matrix_type > inputMatrix_
The input matrix to the constructor.
Definition: Ifpack2_Container_decl.hpp:250
virtual void weightedApplyMV(const mv_type &X, mv_type &Y, vector_type &W) const
Wrapper for weightedApply with MultiVector.
Definition: Ifpack2_Container_def.hpp:150
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:511
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:501
typename mv_type::dual_view_type::t_host HostView
Definition: Ifpack2_Container_decl.hpp:106
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:706
Teuchos::RCP< vector_type > Diag_
Diagonal elements.
Definition: Ifpack2_Container_decl.hpp:267
virtual ~ContainerImpl()
Destructor.
Definition: Ifpack2_Container_def.hpp:493
bool IsParallel_
Whether the problem is distributed across multiple MPI processes.
Definition: Ifpack2_Container_decl.hpp:269
static std::string getName()
Definition: Ifpack2_Container_def.hpp:532
bool IsComputed_
If true, the container has been successfully computed.
Definition: Ifpack2_Container_decl.hpp:289
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
Compute Y := alpha * diag(D) * M^{-1} (diag(D) * X) + beta*Y.
virtual void compute()=0
Extract the local diagonal blocks and prepare the solver.
Interface for creating and solving a set of local linear problems.
Definition: Ifpack2_Container_decl.hpp:79
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:85
virtual void extract()=0
Extract the submatrices identified by the local indices set by the constructor. BlockMatrix may be an...
virtual void setParameters(const Teuchos::ParameterList &List)=0
Set parameters, if any.
Structure for read-only views of general matrix rows.
Definition: Ifpack2_Container_decl.hpp:297
virtual std::ostream & print(std::ostream &os) const =0
Print basic information about the container to os.
static std::string getName()
Definition: Ifpack2_Container_def.hpp:157
GO NumGlobalRows_
Number of global rows in input matrix.
Definition: Ifpack2_Container_decl.hpp:273
GO NumGlobalNonzeros_
Number of nonzeros in input matrix.
Definition: Ifpack2_Container_decl.hpp:275
LO NumLocalRows_
Number of local rows in input matrix.
Definition: Ifpack2_Container_decl.hpp:271
Teuchos::Array< LO > blockOffsets_
Starting index in blockRows_ of local row indices for each block.
Definition: Ifpack2_Container_decl.hpp:265
virtual std::ostream & print(std::ostream &os) const =0
Print basic information about the container to os.
static std::string name()
virtual void solveBlock(ConstHostSubviewLocal X, HostSubviewLocal Y, int blockIndex, Teuchos::ETransp mode, const LSC alpha, const LSC beta) const
Definition: Ifpack2_Container_def.hpp:539
Teuchos::RCP< const crs_matrix_type > inputCrsMatrix_
The input matrix, dynamic cast to CrsMatrix. May be null.
Definition: Ifpack2_Container_decl.hpp:253
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:163