Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_TriDiContainer_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_TRIDICONTAINER_DECL_HPP
10#define IFPACK2_TRIDICONTAINER_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_SerialDenseVector.hpp"
20#include "Teuchos_SerialTriDiMatrix.hpp"
21#include <type_traits>
22#include <string>
23
24namespace Ifpack2 {
25
70
71template <class MatrixType, class LocalScalarType>
73 : public ContainerImpl<MatrixType, LocalScalarType> {
75
76 private:
82 using matrix_type = MatrixType;
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 using LSC = LocalScalarType;
94
95 using typename Container<MatrixType>::mv_type;
96 using typename Container<MatrixType>::map_type;
97 using typename Container<MatrixType>::vector_type;
98 using typename Container<MatrixType>::import_type;
99
100 using typename Container<MatrixType>::ISC;
101 using typename ContainerImpl<MatrixType, LocalScalarType>::LISC;
102 using typename Container<MatrixType>::HostView;
103 using local_mv_type = Tpetra::MultiVector<LSC, LO, GO, NO>;
104 using HostViewLocal = typename Kokkos::View<LSC**, Kokkos::HostSpace>;
105 using typename ContainerImpl<MatrixType, LocalScalarType>::HostSubviewLocal;
106 using typename ContainerImpl<MatrixType, LocalScalarType>::ConstHostSubviewLocal;
107 using typename ContainerImpl<MatrixType, LocalScalarType>::block_crs_matrix_type;
108 static_assert(std::is_same<MatrixType, Tpetra::RowMatrix<SC, LO, GO, NO>>::value,
109 "Ifpack2::TriDiContainer: MatrixType must be a Tpetra::RowMatrix specialization.");
110
119 using typename Container<MatrixType>::row_matrix_type;
121 public:
123
124
135 TriDiContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
136 const Teuchos::Array<Teuchos::Array<LO>>& partitions,
137 const Teuchos::RCP<const import_type>& importer,
138 bool pointIndexed);
139
141
143 virtual ~TriDiContainer();
144
146
148
150
152
154 virtual void initialize();
155
157 virtual void compute();
158
159 void clearBlocks();
160
161 void solveBlock(ConstHostSubviewLocal X,
162 HostSubviewLocal Y,
163 int blockIndex,
164 Teuchos::ETransp mode,
165 LSC alpha,
166 LSC beta) const;
167
169
171
175 virtual std::ostream& print(std::ostream& os) const;
176
178
180
182 virtual std::string description() const;
183
185 virtual void
186 describe(Teuchos::FancyOStream& out,
187 const Teuchos::EVerbosityLevel verbLevel =
188 Teuchos::Describable::verbLevel_default) const;
189
191
193 static std::string getName();
194
195 private:
197 void extract();
198
202 void factor();
203
205 std::vector<Teuchos::SerialTriDiMatrix<int, LSC>> diagBlocks_;
206
208 mutable Teuchos::Array<int> ipiv_;
209
211 Teuchos::Array<LSC> scalars_;
212
214 Teuchos::Array<int> scalarOffsets_;
215};
216
217} // namespace Ifpack2
218
219#endif // IFPACK2_TRIDICONTAINER_DECL_HPP
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
Store and solve a local TriDi linear problem.
Definition Ifpack2_TriDiContainer_decl.hpp:73
virtual std::string description() const
A one-line description of this object.
Definition Ifpack2_TriDiContainer_def.hpp:338
virtual void initialize()
Do all set-up operations that only require matrix structure.
Definition Ifpack2_TriDiContainer_def.hpp:66
virtual void compute()
Extract the local diagonal block and prepare the solver.
Definition Ifpack2_TriDiContainer_def.hpp:76
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_TriDiContainer_def.hpp:357
virtual ~TriDiContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition Ifpack2_TriDiContainer_def.hpp:63
virtual std::ostream & print(std::ostream &os) const
Print information about this object to the given output stream.
Definition Ifpack2_TriDiContainer_def.hpp:330
static std::string getName()
Get the name of this container type for Details::constructContainer()
Definition Ifpack2_TriDiContainer_def.hpp:371
void solveBlock(ConstHostSubviewLocal X, HostSubviewLocal Y, int blockIndex, Teuchos::ETransp mode, LSC alpha, LSC beta) const
Definition Ifpack2_TriDiContainer_def.hpp:226
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40