43 #ifndef IFPACK2_DENSECONTAINER_DECL_HPP
44 #define IFPACK2_DENSECONTAINER_DECL_HPP
49 #include "Ifpack2_Container.hpp"
50 #include "Tpetra_MultiVector.hpp"
51 #include "Tpetra_Map.hpp"
52 #include "Tpetra_RowMatrix.hpp"
54 #include "Tpetra_BlockCrsMatrix_decl.hpp"
102 template<
class MatrixType,
class LocalScalarType>
115 using matrix_type = MatrixType;
117 using LSC = LocalScalarType;
120 using typename Container<MatrixType>::SC;
122 using typename Container<MatrixType>::LO;
124 using typename Container<MatrixType>::GO;
126 using typename Container<MatrixType>::NO;
128 using typename Container<MatrixType>::mv_type;
129 using local_mv_type = Tpetra::MultiVector<LSC, LO, GO, NO>;
130 using typename Container<MatrixType>::map_type;
131 using typename Container<MatrixType>::vector_type;
132 using typename Container<MatrixType>::import_type;
135 using typename ContainerImpl<MatrixType, LocalScalarType>::LISC;
137 using typename ContainerImpl<MatrixType, LocalScalarType>::HostViewLocal;
138 using typename ContainerImpl<MatrixType, LocalScalarType>::HostSubviewLocal;
140 static_assert(std::is_same<MatrixType, Tpetra::RowMatrix<SC, LO, GO, NO>>::value,
141 "Ifpack2::DenseContainer: Please use MatrixType = Tpetra::RowMatrix.");
151 using typename Container<MatrixType>::row_matrix_type;
153 using block_crs_matrix_type = Tpetra::BlockCrsMatrix<SC, LO, GO, NO>;
203 virtual std::ostream&
print (std::ostream& os)
const;
241 solveBlock(HostSubviewLocal X,
246 const LSC beta)
const;
249 std::vector<Teuchos::SerialDenseMatrix<int, LSC>> diagBlocks_;
255 bool hasBlockCrsMatrix_;
266 #endif // IFPACK2_DENSECONTAINER_DECL_HPP
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:342
virtual void compute()
Extract the local diagonal block and prepare the solver.
Definition: Ifpack2_DenseContainer_def.hpp:114
static std::string getName()
Get the name of this container type for Details::constructContainer()
Definition: Ifpack2_DenseContainer_def.hpp:432
DenseContainer(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, const Teuchos::RCP< const import_type > &importer, bool pointIndexed)
Constructor.
Definition: Ifpack2_DenseContainer_def.hpp:61
typename Kokkos::Details::ArithTraits< SC >::val_type ISC
Internal representation of Scalar in Kokkos::View.
Definition: Ifpack2_Container_decl.hpp:135
Store and solve a local dense linear problem.
Definition: Ifpack2_DenseContainer_decl.hpp:103
virtual std::ostream & print(std::ostream &os) const
Print information about this object to the given output stream.
Definition: Ifpack2_DenseContainer_def.hpp:372
typename mv_type::dual_view_type::t_host HostView
Definition: Ifpack2_Container_decl.hpp:139
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to the given FancyOStream.
Definition: Ifpack2_DenseContainer_def.hpp:406
virtual ~DenseContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition: Ifpack2_DenseContainer_def.hpp:93
static const EVerbosityLevel verbLevel_default
virtual std::string description() const
A one-line description of this object.
Definition: Ifpack2_DenseContainer_def.hpp:383
virtual void initialize()
Do all set-up operations that only require matrix structure.
Definition: Ifpack2_DenseContainer_def.hpp:98