Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_Chebyshev_decl.hpp
Go to the documentation of this file.
1 /*
2 @HEADER
3 // ***********************************************************************
4 //
5 // Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
6 // Copyright (2009) Sandia Corporation
7 //
8 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9 // license for use of this work by or on behalf of the U.S. Government.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ***********************************************************************
41 //@HEADER
42 */
43 
50 
51 #ifndef IFPACK2_CHEBYSHEV_DECL_HPP
52 #define IFPACK2_CHEBYSHEV_DECL_HPP
53 
56 
57 // FIXME (mfh 20 Nov 2013) We really shouldn't have to include this.
58 // We should handle the implementation by pointer instead of by value.
59 #include "Ifpack2_Details_Chebyshev.hpp"
60 
61 // We only need the declaration here, and only for the method
62 // getCrsMatrix. I would very much prefer that this method not exist.
63 // Furthermore, it is both unsafe (MatrixType need not be CrsMatrix)
64 // and completely redundant (just call getMatrix() and do the
65 // dynamic_cast yourself).
66 #include "Tpetra_CrsMatrix_decl.hpp"
67 
68 #include <type_traits>
69 
70 namespace Ifpack2 {
71 
198 template<class MatrixType>
199 class Chebyshev :
200  virtual public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
201  typename MatrixType::local_ordinal_type,
202  typename MatrixType::global_ordinal_type,
203  typename MatrixType::node_type>,
204  virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
205  typename MatrixType::local_ordinal_type,
206  typename MatrixType::global_ordinal_type,
207  typename MatrixType::node_type> >
208 {
209 public:
211 
212 
214  typedef MatrixType matrix_type;
215 
217  typedef typename MatrixType::scalar_type scalar_type;
218 
220  typedef typename MatrixType::local_ordinal_type local_ordinal_type;
221 
223  typedef typename MatrixType::global_ordinal_type global_ordinal_type;
224 
226  typedef typename MatrixType::node_type::device_type device_type;
227 
229  typedef typename MatrixType::node_type node_type;
230 
233 
238  typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type,
240 
241  static_assert (std::is_same<MatrixType, row_matrix_type>::value,
242  "Ifpack2::Chebyshev: MatrixType must be a Tpetra::RowMatrix "
243  "specialization. Don't use Tpetra::CrsMatrix here.");
244 
246  typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
247 
253  typedef Tpetra::Vector<scalar_type, local_ordinal_type,
255 
257  // \name Constructors and destructors
259 
282 
284  virtual ~Chebyshev ();
285 
287 
289 
426  void setParameters (const Teuchos::ParameterList& params);
427 
434  void initialize();
435 
438  inline bool isInitialized() const {
439  return IsInitialized_;
440  }
441 
474  void compute ();
475 
485  inline bool isComputed() const {
486  return IsComputed_;
487  }
488 
490 
492 
515  virtual void
517 
519 
521 
550  void
551  apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
552  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
556 
559 
562 
564  bool hasTransposeApply() const;
565 
586  void
587  applyMat (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
588  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
589  Teuchos::ETransp mode = Teuchos::NO_TRANS) const;
590 
592 
594 
597 
600 
606  getCrsMatrix() const;
607 
609  double getComputeFlops() const;
610 
612  double getApplyFlops() const;
613 
615  int getNumInitialize() const;
616 
618  int getNumCompute() const;
619 
621  int getNumApply() const;
622 
624  double getInitializeTime() const;
625 
627  double getComputeTime() const;
628 
630  double getApplyTime() const;
631 
633  size_t getNodeSmootherComplexity() const;
634 
636  typename MatrixType::scalar_type getLambdaMaxForApply() const;
637 
639 
641 
643  std::string description() const;
644 
647 
649 
650 private:
651 
654 
656  typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> MV;
657 
660 
662  Chebyshev<MatrixType>& operator= (const Chebyshev<MatrixType>&);
663 
671  void
672  applyImpl (const MV& X,
673  MV& Y,
674  Teuchos::ETransp mode,
675  scalar_type alpha,
676  scalar_type beta) const;
677 
679 
680 
689 
691  bool IsInitialized_;
693  bool IsComputed_;
695  int NumInitialize_;
697  int NumCompute_;
702  mutable int NumApply_;
704  double InitializeTime_;
706  double ComputeTime_;
711  mutable double ApplyTime_;
713  double ComputeFlops_;
718  mutable double ApplyFlops_;
719 
721 }; // class Chebyshev
722 
723 } // namespace Ifpack2
724 
725 #endif // IFPACK2_CHEBYSHEV_DECL_HPP
726 
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
double getApplyFlops() const
The total number of floating-point operations taken by all calls to apply().
Definition: Ifpack2_Chebyshev_def.hpp:211
Teuchos::RCP< const Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > > getCrsMatrix() const
Attempt to return the matrix A as a Tpetra::CrsMatrix.
Definition: Ifpack2_Chebyshev_def.hpp:127
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Chebyshev_decl.hpp:223
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Chebyshev_decl.hpp:217
double getApplyTime() const
The total time spent in all calls to apply().
Definition: Ifpack2_Chebyshev_def.hpp:199
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Chebyshev_def.hpp:81
Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > vector_type
The Tpetra::Vector specialization matching MatrixType.
Definition: Ifpack2_Chebyshev_decl.hpp:254
void setParameters(const Teuchos::ParameterList &params)
Set (or reset) parameters.
Definition: Ifpack2_Chebyshev_def.hpp:93
double getComputeTime() const
The total time spent in all calls to compute().
Definition: Ifpack2_Chebyshev_def.hpp:193
int getNumApply() const
The total number of successful calls to apply().
Definition: Ifpack2_Chebyshev_def.hpp:181
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Chebyshev_decl.hpp:229
MatrixType::node_type::device_type device_type
The Kokkos::Device specialization used by the input MatrixType.
Definition: Ifpack2_Chebyshev_decl.hpp:226
bool isInitialized() const
Definition: Ifpack2_Chebyshev_decl.hpp:438
void compute()
(Re)compute the left scaling, and (if applicable) estimate max and min eigenvalues of D_inv * A...
Definition: Ifpack2_Chebyshev_def.hpp:305
bool hasTransposeApply() const
Whether it&#39;s possible to apply the transpose of this operator.
Definition: Ifpack2_Chebyshev_def.hpp:163
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
The Tpetra::RowMatrix specialization matching MatrixType.
Definition: Ifpack2_Chebyshev_decl.hpp:239
double getComputeFlops() const
The total number of floating-point operations taken by all calls to compute().
Definition: Ifpack2_Chebyshev_def.hpp:205
Diagonally scaled Chebyshev iteration for Tpetra sparse matrices.
Definition: Ifpack2_Chebyshev_decl.hpp:199
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Apply the preconditioner to X, returning the result in Y.
Definition: Ifpack2_Chebyshev_def.hpp:231
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:107
Chebyshev(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition: Ifpack2_Chebyshev_def.hpp:58
bool isComputed() const
Definition: Ifpack2_Chebyshev_decl.hpp:485
std::string description() const
A simple one-line description of this object.
Definition: Ifpack2_Chebyshev_def.hpp:332
void initialize()
Initialize the preconditioner.
Definition: Ifpack2_Chebyshev_def.hpp:290
Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > map_type
The Tpetra::Map specialization matching MatrixType.
Definition: Ifpack2_Chebyshev_decl.hpp:243
Declaration of interface for preconditioners that can change their matrix after construction.
static const EVerbosityLevel verbLevel_default
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the matrix is distributed.
Definition: Ifpack2_Chebyshev_def.hpp:102
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Chebyshev_decl.hpp:220
int getNumInitialize() const
The total number of successful calls to initialize().
Definition: Ifpack2_Chebyshev_def.hpp:169
double getInitializeTime() const
The total time spent in all calls to initialize().
Definition: Ifpack2_Chebyshev_def.hpp:187
MatrixType::scalar_type getLambdaMaxForApply() const
The estimate of the maximum eigenvalue used in the apply().
Definition: Ifpack2_Chebyshev_def.hpp:499
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Chebyshev_decl.hpp:232
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_Chebyshev_def.hpp:216
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to a Teuchos::FancyOStream.
Definition: Ifpack2_Chebyshev_def.hpp:361
Teuchos::RCP< const row_matrix_type > getMatrix() const
The matrix for which this is a preconditioner.
Definition: Ifpack2_Chebyshev_def.hpp:116
Teuchos::RCP< const map_type > getDomainMap() const
The Tpetra::Map representing the domain of this operator.
Definition: Ifpack2_Chebyshev_def.hpp:137
MatrixType matrix_type
The template parameter of this class.
Definition: Ifpack2_Chebyshev_decl.hpp:214
virtual ~Chebyshev()
Destructor.
Definition: Ifpack2_Chebyshev_def.hpp:76
Teuchos::RCP< const map_type > getRangeMap() const
The Tpetra::Map representing the range of this operator.
Definition: Ifpack2_Chebyshev_def.hpp:151
int getNumCompute() const
The total number of successful calls to compute().
Definition: Ifpack2_Chebyshev_def.hpp:175
void applyMat(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS) const
Compute Y = Op(A)*X, where Op(A) is either A, , or .
Definition: Ifpack2_Chebyshev_def.hpp:271