Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_BandedContainer_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9#ifndef IFPACK2_BANDEDCONTAINER_DECL_HPP
10#define IFPACK2_BANDEDCONTAINER_DECL_HPP
11
14
15#include "Ifpack2_Container.hpp"
16#include "Tpetra_MultiVector.hpp"
17#include "Tpetra_Map.hpp"
18#include "Tpetra_RowMatrix.hpp"
19#include "Teuchos_SerialBandDenseMatrix.hpp"
20
21namespace Ifpack2 {
22
67
68template<class MatrixType, class LocalScalarType>
70: public ContainerImpl<MatrixType, LocalScalarType>
71{
73
74private:
81 using matrix_type = MatrixType;
83 using LSC = LocalScalarType;
85
87 using typename Container<MatrixType>::SC;
89 using typename Container<MatrixType>::LO;
91 using typename Container<MatrixType>::GO;
93 using typename Container<MatrixType>::NO;
94
95 using typename Container<MatrixType>::ISC;
96 using typename ContainerImpl<MatrixType, LSC>::LISC;
97
98 using typename Container<MatrixType>::mv_type;
99 using typename Container<MatrixType>::map_type;
100 using local_mv_type = Tpetra::MultiVector<LSC, LO, GO, NO>;
101 using typename Container<MatrixType>::vector_type;
102 using typename Container<MatrixType>::import_type;
103
104 using typename Container<MatrixType>::HostView;
105 using typename ContainerImpl<MatrixType, LSC>::HostSubviewLocal;
106 using typename ContainerImpl<MatrixType, LSC>::ConstHostSubviewLocal;
107 using typename ContainerImpl<MatrixType,LSC>::block_crs_matrix_type;
108 using HostViewLocal = typename local_mv_type::dual_view_type::t_host;
109
110 static_assert(std::is_same<MatrixType,
111 Tpetra::RowMatrix<SC, LO, GO, NO> >::value,
112 "Ifpack2::BandedContainer: Please use MatrixType = Tpetra::RowMatrix.");
113
122 using typename Container<MatrixType>::row_matrix_type;
123
125public:
127
128
139 BandedContainer (const Teuchos::RCP<const row_matrix_type>& matrix,
140 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
141 const Teuchos::RCP<const import_type>&,
142 bool pointIndexed);
143
145
147 virtual ~BandedContainer ();
148
150
152
154 virtual void setParameters (const Teuchos::ParameterList& List) override;
155
157
159
161 virtual void initialize () override;
162
164 virtual void compute () override;
165
168 void computeBandwidth();
169
170 void clearBlocks() override;
171
173
175
179 virtual std::ostream& print (std::ostream& os) const override;
180
182
184
186 virtual std::string description () const override;
187
189 virtual void
190 describe (Teuchos::FancyOStream &out,
191 const Teuchos::EVerbosityLevel verbLevel =
192 Teuchos::Describable::verbLevel_default) const override;
193
195
197 static std::string getName();
198
199private:
200
202 void extract() override;
203
207 void factor ();
208
217 void
218 solveBlock(ConstHostSubviewLocal X,
219 HostSubviewLocal Y,
220 int blockIndex,
221 Teuchos::ETransp mode,
222 const LSC alpha,
223 const LSC beta) const override;
224
226 std::vector<Teuchos::SerialBandDenseMatrix<int, LSC> > diagBlocks_;
227
229 Teuchos::Array<int> ipiv_;
230
231 Teuchos::Array<LO> kl_; //< number of subdiagonals
232 Teuchos::Array<LO> ku_; //< number of superdiagonals
233
235 Teuchos::Array<LSC> scalars_;
236
238 Teuchos::Array<LO> scalarOffsets_;
239};
240
241} // namespace Ifpack2
242
243#endif // IFPACK2_BANDEDCONTAINER_DECL_HPP
BandedContainer(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, const Teuchos::RCP< const import_type > &, bool pointIndexed)
Constructor.
Definition Ifpack2_BandedContainer_def.hpp:29
virtual void compute() override
Initialize and compute each block.
Definition Ifpack2_BandedContainer_def.hpp:209
virtual void initialize() override
Do all set-up operations that only require matrix structure.
Definition Ifpack2_BandedContainer_def.hpp:175
virtual std::ostream & print(std::ostream &os) const override
Print information about this object to the given output stream.
Definition Ifpack2_BandedContainer_def.hpp:497
virtual ~BandedContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition Ifpack2_BandedContainer_def.hpp:46
virtual void setParameters(const Teuchos::ParameterList &List) override
Set all necessary parameters (super- and sub-diagonal bandwidth).
Definition Ifpack2_BandedContainer_def.hpp:50
static std::string getName()
Get the name of this container type for Details::constructContainer().
Definition Ifpack2_BandedContainer_def.hpp:549
virtual std::string description() const override
A one-line description of this object.
Definition Ifpack2_BandedContainer_def.hpp:508
void computeBandwidth()
Definition Ifpack2_BandedContainer_def.hpp:68
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Print the object with some verbosity level to the given FancyOStream.
Definition Ifpack2_BandedContainer_def.hpp:530
virtual void extract()=0
Extract the submatrices identified by the local indices set by the constructor. BlockMatrix may be an...
typename Kokkos::ArithTraits< SC >::val_type ISC
Internal representation of Scalar in Kokkos::View.
Definition Ifpack2_Container_decl.hpp:102
typename mv_type::dual_view_type::t_host HostView
Definition Ifpack2_Container_decl.hpp:106
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41