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>
21#include <Teuchos_ParameterList.hpp>
24#ifndef DOXYGEN_SHOULD_SKIP_THIS
78template <
class MatrixType>
81 using scalar_type =
typename MatrixType::scalar_type;
82 using local_ordinal_type =
typename MatrixType::local_ordinal_type;
83 using global_ordinal_type =
typename MatrixType::global_ordinal_type;
84 using node_type =
typename MatrixType::node_type;
85 using SC = scalar_type;
86 using LO = local_ordinal_type;
87 using GO = global_ordinal_type;
89 using import_type = Tpetra::Import<LO, GO, NO>;
90 using mv_type = Tpetra::MultiVector<SC, LO, GO, NO>;
91 using vector_type = Tpetra::Vector<SC, LO, GO, NO>;
92 using map_type = Tpetra::Map<LO, GO, NO>;
93 using crs_matrix_type = Tpetra::CrsMatrix<SC, LO, GO, NO>;
94 using block_crs_matrix_type = Tpetra::BlockCrsMatrix<SC, LO, GO, NO>;
95 using row_matrix_type = Tpetra::RowMatrix<SC, LO, GO, NO>;
97 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
98 "Ifpack2::Container: Please use MatrixType = Tpetra::RowMatrix.");
101#if KOKKOS_VERSION >= 40799
102 using ISC =
typename KokkosKernels::ArithTraits<SC>::val_type;
104 using ISC =
typename Kokkos::ArithTraits<SC>::val_type;
109 using HostView =
typename mv_type::dual_view_type::t_host;
110 using ConstHostView =
typename HostView::const_type;
123 Container(
const Teuchos::RCP<const row_matrix_type>& matrix,
124 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
145 Teuchos::ArrayView<const LO>
getBlockRows(
int blockIndex)
const;
158 void setBlockSizes(
const Teuchos::Array<Teuchos::Array<LO> >& partitions);
160 void getMatDiag()
const;
178 void DoJacobi(ConstHostView X,
HostView Y, SC dampingFactor)
const;
179 void DoOverlappingJacobi(ConstHostView X,
HostView Y, ConstHostView W, SC dampingFactor,
bool nonsymScaling)
const;
180 void DoGaussSeidel(ConstHostView X,
HostView Y,
HostView Y2, SC dampingFactor)
const;
181 void DoSGS(ConstHostView X,
HostView Y,
HostView Y2, SC dampingFactor)
const;
202 Teuchos::ETransp mode = Teuchos::NO_TRANS,
203 SC alpha = Teuchos::ScalarTraits<SC>::one(),
204 SC beta = Teuchos::ScalarTraits<SC>::zero())
const = 0;
212 Teuchos::ETransp mode = Teuchos::NO_TRANS,
213 SC alpha = Teuchos::ScalarTraits<SC>::one(),
214 SC beta = Teuchos::ScalarTraits<SC>::zero())
const = 0;
230 virtual void applyMV(
const mv_type& X, mv_type& Y)
const;
235 vector_type& W)
const;
237 virtual void clearBlocks();
240 virtual std::ostream&
print(std::ostream& os)
const = 0;
249 SC dampingFactor, LO i)
const;
269 mutable Teuchos::RCP<vector_type>
Diag_;
297template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
298struct StridedRowView;
310template <
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#if KOKKOS_VERSION >= 40799
332 using LISC =
typename KokkosKernels::ArithTraits<LSC>::val_type;
334 using LISC =
typename Kokkos::ArithTraits<LSC>::val_type;
337 using local_mv_type = Tpetra::MultiVector<LSC, LO, GO, NO>;
340 using typename Container<MatrixType>::ConstHostView;
341 using HostViewLocal =
typename local_mv_type::dual_view_type::t_host;
342 using HostSubviewLocal = Kokkos::View<LISC**, Kokkos::LayoutStride, typename HostViewLocal::memory_space>;
343 using ConstHostSubviewLocal = Kokkos::View<const LISC**, Kokkos::LayoutStride, typename HostViewLocal::memory_space>;
345 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
346 "Ifpack2::Container: Please use MatrixType = Tpetra::RowMatrix.");
368 ContainerImpl(
const Teuchos::RCP<const row_matrix_type>& matrix,
369 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
396 virtual void setParameters(
const Teuchos::ParameterList& List);
411 apply(ConstHostView X,
414 Teuchos::ETransp mode = Teuchos::NO_TRANS,
415 SC alpha = Teuchos::ScalarTraits<SC>::one(),
416 SC beta = Teuchos::ScalarTraits<SC>::zero())
const;
424 Teuchos::ETransp mode = Teuchos::NO_TRANS,
425 SC alpha = Teuchos::ScalarTraits<SC>::one(),
426 SC beta = Teuchos::ScalarTraits<SC>::zero())
const;
442 void applyMV(
const mv_type& X, mv_type& Y)
const;
447 vector_type& W)
const;
449 virtual void clearBlocks();
452 virtual std::ostream&
print(std::ostream& os)
const = 0;
461 SC dampingFactor, LO i)
const;
472 Teuchos::ETransp mode,
474 const LSC beta)
const;
478 mutable HostViewLocal Y_local_;
505template <
typename Scalar,
typename LocalOrdinal,
typename GlobalOrdinal,
typename Node>
508 using LO = LocalOrdinal;
510 using block_crs_matrix_type = Tpetra::BlockCrsMatrix<SC, LO, GlobalOrdinal, Node>;
512 using h_inds_type =
typename block_crs_matrix_type::local_inds_host_view_type;
513 using h_vals_type =
typename block_crs_matrix_type::values_host_view_type;
523 SC val(
size_t i)
const;
524 LO ind(
size_t i)
const;
534 Teuchos::Array<SC> valsCopy;
535 Teuchos::Array<LO> indsCopy;
542template <
class MatrixType>
551template <
class MatrixType>
554 static std::string name() {
555 return std::string(
"Ifpack2::Container<") +
556 TypeNameTraits<MatrixType>::name() +
">";
559 static std::string concreteName(const ::Ifpack2::Container<MatrixType>&) {
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS)....
Definition Ifpack2_Container_decl.hpp:311
std::vector< HostSubviewLocal > X_localBlocks_
Views for holding pieces of X corresponding to each block.
Definition Ifpack2_Container_decl.hpp:493
HostViewLocal X_local_
Scratch vectors used in apply().
Definition Ifpack2_Container_decl.hpp:477
virtual std::ostream & print(std::ostream &os) const =0
Print basic information about the container to os.
virtual void setParameters(const Teuchos::ParameterList &List)
Set parameters, if any.
Definition Ifpack2_Container_def.hpp:443
Details::StridedRowView< SC, LO, GO, NO > getInputRowView(LO row) const
View a row of the input matrix.
Definition Ifpack2_Container_def.hpp:803
virtual void extract()=0
Extract the submatrices identified by the local indices set by the constructor. BlockMatrix may be an...
std::vector< HostSubviewLocal > Y_localBlocks_
Views for holding pieces of Y corresponding to each block.
Definition Ifpack2_Container_decl.hpp:496
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:648
virtual void compute()=0
Extract the local diagonal blocks and prepare the solver.
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:456
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:539
LocalScalarType LSC
The internal representation of LocalScalarType in Kokkos::View.
Definition Ifpack2_Container_decl.hpp:330
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:231
virtual void solveBlock(ConstHostSubviewLocal X, HostSubviewLocal Y, int blockIndex, Teuchos::ETransp mode, const LSC alpha, const LSC beta) const
Definition Ifpack2_Container_def.hpp:481
HostViewLocal weightedApplyScratch_
Definition Ifpack2_Container_decl.hpp:490
static std::string getName()
Definition Ifpack2_Container_def.hpp:475
LO translateRowToCol(LO row)
Definition Ifpack2_Container_def.hpp:493
virtual void initialize()=0
Do all set-up operations that only require matrix structure.
virtual ~ContainerImpl()
Destructor.
Definition Ifpack2_Container_def.hpp:439
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:464
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:447
Interface for creating and solving a set of local linear problems.
Definition Ifpack2_Container_decl.hpp:79
virtual ~Container()
Destructor.
Definition Ifpack2_Container_def.hpp:80
bool hasBlockCrs_
Whether the input matrix is a BlockCRS matrix.
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.
Teuchos::RCP< const crs_matrix_type > inputCrsMatrix_
The input matrix, dynamic cast to CrsMatrix. May be null.
Definition Ifpack2_Container_decl.hpp:255
bool IsParallel_
Whether the problem is distributed across multiple MPI processes.
Definition Ifpack2_Container_decl.hpp:271
int numBlocks_
The number of blocks (partitions) in the container.
Definition Ifpack2_Container_decl.hpp:261
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:84
static std::string getName()
Definition Ifpack2_Container_def.hpp:146
typename Kokkos::ArithTraits< SC >::val_type ISC
Internal representation of Scalar in Kokkos::View.
Definition Ifpack2_Container_decl.hpp:104
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).
Teuchos::RCP< vector_type > Diag_
Diagonal elements.
Definition Ifpack2_Container_decl.hpp:269
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.
int bcrsBlockSize_
If hasBlockCrs_, the number of DOFs per vertex. Otherwise 1.
Definition Ifpack2_Container_decl.hpp:281
Teuchos::RCP< const block_crs_matrix_type > inputBlockMatrix_
The input matrix, dynamic cast to BlockCrsMatrix. May be null.
Definition Ifpack2_Container_decl.hpp:258
virtual void setParameters(const Teuchos::ParameterList &List)=0
Set parameters, if any.
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:151
virtual void applyMV(const mv_type &X, mv_type &Y) const
Wrapper for apply with MultiVector.
Definition Ifpack2_Container_def.hpp:134
LO NumLocalRows_
Number of local rows in input matrix.
Definition Ifpack2_Container_decl.hpp:273
bool isInitialized() const
Whether the container has been successfully initialized.
Definition Ifpack2_Container_def.hpp:123
void setBlockSizes(const Teuchos::Array< Teuchos::Array< LO > > &partitions)
Initialize arrays with information about block sizes.
Definition Ifpack2_Container_def.hpp:89
LO scalarsPerRow_
Definition Ifpack2_Container_decl.hpp:286
Teuchos::Array< LO > blockSizes_
Number of rows in each block.
Definition Ifpack2_Container_decl.hpp:265
GO NumGlobalNonzeros_
Number of nonzeros in input matrix.
Definition Ifpack2_Container_decl.hpp:277
virtual void initialize()=0
Do all set-up operations that only require matrix structure.
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:140
Teuchos::RCP< const row_matrix_type > inputMatrix_
The input matrix to the constructor.
Definition Ifpack2_Container_decl.hpp:252
bool pointIndexed_
(If hasBlockCrs_) Whether the blocks are described using sub-block row indices instead of full block ...
Definition Ifpack2_Container_decl.hpp:283
bool IsInitialized_
If true, the container has been successfully initialized.
Definition Ifpack2_Container_decl.hpp:289
Teuchos::Array< LO > blockRows_
Local indices of the rows of the input matrix that belong to this block.
Definition Ifpack2_Container_decl.hpp:263
GO NumGlobalRows_
Number of global rows in input matrix.
Definition Ifpack2_Container_decl.hpp:275
bool isComputed() const
Whether the container has been successfully computed.
Definition Ifpack2_Container_def.hpp:128
typename mv_type::dual_view_type::t_host HostView
Definition Ifpack2_Container_decl.hpp:109
Teuchos::Array< LO > blockOffsets_
Starting index in blockRows_ of local row indices for each block.
Definition Ifpack2_Container_decl.hpp:267
bool IsComputed_
If true, the container has been successfully computed.
Definition Ifpack2_Container_decl.hpp:291
virtual void compute()=0
Extract the local diagonal blocks and prepare the solver.
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40
Structure for read-only views of general matrix rows.
Definition Ifpack2_Container_decl.hpp:506