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> {
72
73 private:
80 using matrix_type = MatrixType;
82 using LSC = LocalScalarType;
84
86 using typename Container<MatrixType>::SC;
88 using typename Container<MatrixType>::LO;
90 using typename Container<MatrixType>::GO;
92 using typename Container<MatrixType>::NO;
93
94 using typename Container<MatrixType>::ISC;
95 using typename ContainerImpl<MatrixType, LSC>::LISC;
96
97 using typename Container<MatrixType>::mv_type;
98 using typename Container<MatrixType>::map_type;
99 using local_mv_type = Tpetra::MultiVector<LSC, LO, GO, NO>;
100 using typename Container<MatrixType>::vector_type;
101 using typename Container<MatrixType>::import_type;
102
103 using typename Container<MatrixType>::HostView;
104 using typename ContainerImpl<MatrixType, LSC>::HostSubviewLocal;
105 using typename ContainerImpl<MatrixType, LSC>::ConstHostSubviewLocal;
106 using typename ContainerImpl<MatrixType, LSC>::block_crs_matrix_type;
107 using HostViewLocal = typename local_mv_type::dual_view_type::t_host;
108
109 static_assert(std::is_same<MatrixType,
110 Tpetra::RowMatrix<SC, LO, GO, NO> >::value,
111 "Ifpack2::BandedContainer: Please use MatrixType = Tpetra::RowMatrix.");
112
121 using typename Container<MatrixType>::row_matrix_type;
122
124 public:
126
127
138 BandedContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
139 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
140 const Teuchos::RCP<const import_type>&,
141 bool pointIndexed);
142
144
146 virtual ~BandedContainer();
147
149
151
153 virtual void setParameters(const Teuchos::ParameterList& List) override;
154
156
158
160 virtual void initialize() override;
161
163 virtual void compute() override;
164
167 void computeBandwidth();
168
169 void clearBlocks() override;
170
172
174
178 virtual std::ostream& print(std::ostream& os) const override;
179
181
183
185 virtual std::string description() const override;
186
188 virtual void
189 describe(Teuchos::FancyOStream& out,
190 const Teuchos::EVerbosityLevel verbLevel =
191 Teuchos::Describable::verbLevel_default) const override;
192
194
196 static std::string getName();
197
198 private:
200 void extract() override;
201
205 void factor();
206
215 void
216 solveBlock(ConstHostSubviewLocal X,
217 HostSubviewLocal Y,
218 int blockIndex,
219 Teuchos::ETransp mode,
220 const LSC alpha,
221 const LSC beta) const override;
222
224 std::vector<Teuchos::SerialBandDenseMatrix<int, LSC> > diagBlocks_;
225
227 Teuchos::Array<int> ipiv_;
228
229 Teuchos::Array<LO> kl_; //< number of subdiagonals
230 Teuchos::Array<LO> ku_; //< number of superdiagonals
231
233 Teuchos::Array<LSC> scalars_;
234
236 Teuchos::Array<LO> scalarOffsets_;
237};
238
239} // namespace Ifpack2
240
241#endif // IFPACK2_BANDEDCONTAINER_DECL_HPP
Store and solve a local Banded linear problem.
Definition Ifpack2_BandedContainer_decl.hpp:70
virtual void compute() override
Initialize and compute each block.
Definition Ifpack2_BandedContainer_def.hpp:186
virtual void initialize() override
Do all set-up operations that only require matrix structure.
Definition Ifpack2_BandedContainer_def.hpp:156
virtual std::ostream & print(std::ostream &os) const override
Print information about this object to the given output stream.
Definition Ifpack2_BandedContainer_def.hpp:455
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:500
virtual std::string description() const override
A one-line description of this object.
Definition Ifpack2_BandedContainer_def.hpp:465
void computeBandwidth()
Definition Ifpack2_BandedContainer_def.hpp:65
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:483
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS)....
Definition Ifpack2_Container_decl.hpp:311
typename Kokkos::ArithTraits< SC >::val_type ISC
Internal representation of Scalar in Kokkos::View.
Definition Ifpack2_Container_decl.hpp:104
typename mv_type::dual_view_type::t_host HostView
Definition Ifpack2_Container_decl.hpp:109
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40