Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_Details_TriDiSolver_decl.hpp
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 
10 #ifndef IFPACK2_DETAILS_TRIDISOLVER_DECL_HPP
11 #define IFPACK2_DETAILS_TRIDISOLVER_DECL_HPP
12 
15 
16 #include "Ifpack2_ConfigDefs.hpp"
19 #include "Ifpack2_Details_LapackSupportsScalar.hpp"
20 #include "Tpetra_RowMatrix.hpp"
21 #include "Tpetra_Import_fwd.hpp"
22 #include "Tpetra_Export_fwd.hpp"
24 #include <type_traits>
25 
26 namespace Ifpack2 {
27 namespace Details {
28 
41 template<class MatrixType,
42  const bool stub = ! LapackSupportsScalar<typename MatrixType::scalar_type>::value>
43 class TriDiSolver :
44  public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
45  typename MatrixType::local_ordinal_type,
46  typename MatrixType::global_ordinal_type,
47  typename MatrixType::node_type>,
48  virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
49  typename MatrixType::local_ordinal_type,
50  typename MatrixType::global_ordinal_type,
51  typename MatrixType::node_type> >
52 {};
53 
55 template<class MatrixType>
56 class TriDiSolver<MatrixType, false> :
57  public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
58  typename MatrixType::local_ordinal_type,
59  typename MatrixType::global_ordinal_type,
60  typename MatrixType::node_type>,
61  virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
62  typename MatrixType::local_ordinal_type,
63  typename MatrixType::global_ordinal_type,
64  typename MatrixType::node_type> >
65 {
66 public:
68 
69 
73  typedef MatrixType matrix_type;
74 
76  typedef typename MatrixType::scalar_type scalar_type;
77 
79  typedef typename MatrixType::local_ordinal_type local_ordinal_type;
80 
82  typedef typename MatrixType::global_ordinal_type global_ordinal_type;
83 
85  typedef typename MatrixType::node_type node_type;
86 
89 
91  typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> row_matrix_type;
92 
93  static_assert(std::is_same<MatrixType, row_matrix_type>::value, "Ifpack2::Details::TriDiSolver: The template parameter MatrixType must be a Tpetra::RowMatrix specialization. Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore. The constructor can take either a RowMatrix or a CrsMatrix just fine.");
94 
95  typedef typename row_matrix_type::nonconst_global_inds_host_view_type nonconst_global_inds_host_view_type;
96  typedef typename row_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
97  typedef typename row_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
98 
100  typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
101 
103 
105 
110 
112  virtual ~TriDiSolver ();
113 
115 
117 
118 
124 
130 
136  void
137  apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
138  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
142 
144 
146  void setParameters (const Teuchos::ParameterList& params);
147 
157  void initialize ();
158 
160  bool isInitialized () const;
161 
171  void compute ();
172 
174  bool isComputed () const;
175 
178 
181 
183  int getNumInitialize () const;
184 
186  int getNumCompute () const;
187 
189  int getNumApply() const;
190 
192  double getInitializeTime() const;
193 
195  double getComputeTime() const;
196 
198  double getApplyTime() const;
199 
201 
203 
205  std::string description () const;
206 
208  void
209  describe (Teuchos::FancyOStream &out,
210  const Teuchos::EVerbosityLevel verbLevel =
213 private:
214 
216  void
217  describeLocal (Teuchos::FancyOStream& out,
218  const Teuchos::EVerbosityLevel verbLevel) const;
219 
221  void reset ();
222 
230  static void
231  extract (Teuchos::SerialTriDiMatrix<int, scalar_type>& A_local_tridi,
232  const row_matrix_type& A_local);
233 
243  static void
245  const Teuchos::ArrayView<int>& ipiv);
246 
248  typedef Tpetra::MultiVector<scalar_type, local_ordinal_type,
250 
252  typedef Tpetra::Import<local_ordinal_type,
253  global_ordinal_type, node_type> import_type;
254 
256  typedef Tpetra::Export<local_ordinal_type,
257  global_ordinal_type, node_type> export_type;
258 
261 
270  void
271  applyImpl (const MV& X,
272  MV& Y,
273  const Teuchos::ETransp mode,
274  const scalar_type alpha,
275  const scalar_type beta) const;
276 
279 
282 
285 
287  Teuchos::Array<int> ipiv_;
288 
290  double initializeTime_;
291 
293  double computeTime_;
294 
296  mutable double applyTime_;
297 
299  int numInitialize_;
300 
302  int numCompute_;
303 
305  mutable int numApply_;
306 
308  bool isInitialized_;
309 
311  bool isComputed_;
312 };
313 
314 
316 template<class MatrixType>
317 class TriDiSolver<MatrixType, true> :
318  public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
319  typename MatrixType::local_ordinal_type,
320  typename MatrixType::global_ordinal_type,
321  typename MatrixType::node_type>,
322  virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
323  typename MatrixType::local_ordinal_type,
324  typename MatrixType::global_ordinal_type,
325  typename MatrixType::node_type> >
326 {
327 public:
329 
330 
334  typedef MatrixType matrix_type;
335 
337  typedef typename MatrixType::scalar_type scalar_type;
338 
340  typedef typename MatrixType::local_ordinal_type local_ordinal_type;
341 
343  typedef typename MatrixType::global_ordinal_type global_ordinal_type;
344 
346  typedef typename MatrixType::node_type node_type;
347 
350 
352  typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> row_matrix_type;
353 
354  static_assert(std::is_same<MatrixType, row_matrix_type>::value, "Ifpack2::Details::TriDiSolver: The template parameter MatrixType must be a Tpetra::RowMatrix specialization. Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore. The constructor can take either a RowMatrix or a CrsMatrix just fine.");
355 
357  typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
358 
360 
362 
367 
369  virtual ~TriDiSolver ();
370 
372 
374 
380 
386 
392  void
393  apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
394  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
396  scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
397  scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
398 
400 
402  void setParameters (const Teuchos::ParameterList& params);
403 
413  void initialize ();
414 
416  bool isInitialized () const;
417 
427  void compute ();
428 
430  bool isComputed () const;
431 
434 
437 
439  int getNumInitialize () const;
440 
442  int getNumCompute () const;
443 
445  int getNumApply() const;
446 
448  double getInitializeTime() const;
449 
451  double getComputeTime() const;
452 
454  double getApplyTime() const;
455 
457 
459 
461  std::string description () const;
462 
464  void
465  describe (Teuchos::FancyOStream &out,
466  const Teuchos::EVerbosityLevel verbLevel =
468 
469  void
470  describeLocal (Teuchos::FancyOStream& out,
471  const Teuchos::EVerbosityLevel verbLevel) const;
472 
474 private:
476  typedef Tpetra::MultiVector<scalar_type, local_ordinal_type,
477  global_ordinal_type, node_type> MV;
478 };
479 
480 } // namespace Details
481 } // namespace Ifpack2
482 
483 #endif // IFPACK2_DETAILS_TRIDISOLVER_DECL_HPP
virtual Teuchos::RCP< const Tpetra::Map< MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > getDomainMap() const =0
The domain Map of this operator.
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:60
virtual Teuchos::RCP< const Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > getMatrix() const =0
The input matrix given to the constructor.
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the absolute value (magnitude) of a scalar_type.
Definition: Ifpack2_Details_TriDiSolver_decl.hpp:88
&quot;Preconditioner&quot; that uses LAPACK&#39;s tridi LU.
Definition: Ifpack2_Details_TriDiSolver_decl.hpp:43
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input (global) matrix.
Definition: Ifpack2_Details_TriDiSolver_decl.hpp:82
Teuchos::ScalarTraits< MatrixType::scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Preconditioner.hpp:78
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Specialization of Tpetra::RowMatrix used by this class.
Definition: Ifpack2_Details_TriDiSolver_decl.hpp:91
MatrixType::scalar_type scalar_type
The type of entries in the input (global) matrix.
Definition: Ifpack2_Details_TriDiSolver_decl.hpp:76
virtual Teuchos::RCP< const Tpetra::Map< MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > getRangeMap() const =0
The range Map of this operator.
MatrixType matrix_type
The first template parameter of this class.
Definition: Ifpack2_Details_TriDiSolver_decl.hpp:73
Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > map_type
Specialization of Tpetra::Map used by this class.
Definition: Ifpack2_Details_TriDiSolver_decl.hpp:100
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:74
Declaration of interface for preconditioners that can change their matrix after construction.
static const EVerbosityLevel verbLevel_default
virtual bool isComputed() const =0
True if the preconditioner has been successfully computed, else false.
virtual void setParameters(const Teuchos::ParameterList &List)=0
Set this preconditioner&#39;s parameters.
virtual void apply(const Tpetra::MultiVector< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > &X, Tpetra::MultiVector< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, MatrixType::scalar_typealpha=Teuchos::ScalarTraits< MatrixType::scalar_type >::one(), MatrixType::scalar_typebeta=Teuchos::ScalarTraits< MatrixType::scalar_type >::zero()) const =0
Apply the preconditioner to X, putting the result in Y.
virtual bool isInitialized() const =0
True if the preconditioner has been successfully initialized, else false.
MatrixType::node_type node_type
The Node type of the input (global) matrix.
Definition: Ifpack2_Details_TriDiSolver_decl.hpp:85
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input (global) matrix.
Definition: Ifpack2_Details_TriDiSolver_decl.hpp:79
virtual void setMatrix(const Teuchos::RCP< const Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > &A)=0
Set the new matrix.