Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_AdditiveSchwarz_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 
18 
19 #ifndef IFPACK2_ADDITIVESCHWARZ_DECL_HPP
20 #define IFPACK2_ADDITIVESCHWARZ_DECL_HPP
21 
22 #include "Ifpack2_ConfigDefs.hpp"
26 #include "Ifpack2_OverlappingRowMatrix.hpp"
27 #include "Tpetra_Map.hpp"
28 #include "Tpetra_MultiVector.hpp"
29 #include "Tpetra_RowMatrix.hpp"
30 #include <memory>
31 #include <type_traits>
32 
33 namespace Trilinos {
34 namespace Details {
35 template<class MV, class OP, class NormType>
36 class LinearSolver; // forward declaration
37 } // namespace Details
38 } // namespace Trilinos
39 
40 
41 namespace Ifpack2 {
42 
244 template<class MatrixType,
245  class LocalInverseType =
246  Preconditioner<typename MatrixType::scalar_type,
247  typename MatrixType::local_ordinal_type,
248  typename MatrixType::global_ordinal_type,
249  typename MatrixType::node_type> >
251  virtual public Preconditioner<typename MatrixType::scalar_type,
252  typename MatrixType::local_ordinal_type,
253  typename MatrixType::global_ordinal_type,
254  typename MatrixType::node_type>,
255  virtual public Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
256  typename MatrixType::local_ordinal_type,
257  typename MatrixType::global_ordinal_type,
258  typename MatrixType::node_type> >,
259  virtual public Details::NestedPreconditioner<Preconditioner<typename MatrixType::scalar_type,
260  typename MatrixType::local_ordinal_type,
261  typename MatrixType::global_ordinal_type,
262  typename MatrixType::node_type> >
263 {
264 public:
265  static_assert(std::is_same<LocalInverseType,
266  Preconditioner<typename MatrixType::scalar_type,
267  typename MatrixType::local_ordinal_type,
268  typename MatrixType::global_ordinal_type,
269  typename MatrixType::node_type> >::value, "Ifpack2::AdditiveSchwarz: You are not allowed to use nondefault values for the LocalInverseType template parameter. Please stop specifying this explicitly. The default template parameter is perfectly fine.");
270 
271  static_assert(std::is_same<MatrixType,
272  Tpetra::RowMatrix<typename MatrixType::scalar_type,
273  typename MatrixType::local_ordinal_type,
274  typename MatrixType::global_ordinal_type,
275  typename MatrixType::node_type> >::value, "Ifpack2::AdditiveSchwarz: Please use MatrixType = Tpetra::RowMatrix instead of MatrixType = Tpetra::CrsMatrix. Don't worry, AdditiveSchwarz's constructor can take either type of matrix; it does a dynamic cast if necessary inside. Restricting the set of allowed types here will improve build times and reduce library and executable sizes.");
276 
278 
279 
281  using scalar_type = typename MatrixType::scalar_type;
282 
284  using local_ordinal_type = typename MatrixType::local_ordinal_type;
285 
287  using global_ordinal_type = typename MatrixType::global_ordinal_type;
288 
290  using node_type = typename MatrixType::node_type;
291 
293  using magnitude_type =
295 
297  using row_matrix_type =
298  Tpetra::RowMatrix<scalar_type, local_ordinal_type,
300 
302  using crs_matrix_type =
303  Tpetra::CrsMatrix<scalar_type, local_ordinal_type,
305 
307  // \name Constructors and destructor
309 
314 
325  const int overlapLevel);
326 
328  virtual ~AdditiveSchwarz () = default;
329 
331 
333 
336 
339 
341  virtual void
342  apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &X,
343  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &Y,
347 
349 
351 
386  virtual void
390  node_type> >& innerPrec);
391 
393 
395 
418  virtual void
421 
424 
594  virtual void setParameters (const Teuchos::ParameterList& plist);
595 
619  void
621 
622  bool supportsZeroStartingSolution() { return true; }
623 
624  void setZeroStartingSolution (bool zeroStartingSolution) { ZeroStartingSolution_ = zeroStartingSolution; };
625 
631 
633  virtual void initialize();
634 
636  virtual bool isInitialized() const;
637 
639  virtual void compute();
640 
642  virtual bool isComputed() const;
643 
645  virtual int getNumInitialize() const;
646 
648  virtual int getNumCompute() const;
649 
651  virtual int getNumApply() const;
652 
654  virtual double getInitializeTime() const;
655 
657  virtual double getComputeTime() const;
658 
660  virtual double getApplyTime() const;
661 
663 
664 
666  std::string description() const;
667 
670 
672 
674  virtual std::ostream& print(std::ostream& os) const;
675 
677  virtual int getOverlapLevel() const;
678 
679 
680 private:
682  typedef Tpetra::Map<local_ordinal_type,
684  node_type> map_type;
686  typedef Tpetra::Import<local_ordinal_type,
688  node_type> import_type;
690  typedef Tpetra::MultiVector<scalar_type,
693  node_type> MV;
695  typedef Tpetra::Operator<scalar_type,
698  node_type> OP;
700  typedef Preconditioner<scalar_type,
703  node_type> prec_type;
704 
707 
709  AdditiveSchwarz (const AdditiveSchwarz& RHS);
710 
712  void setup ();
713 
715  void localApply(MV &OverlappingX, MV &OverlappingY) const;
716 
719  bool hasInnerPrecName () const;
720 
721 
723  std::string innerPrecName () const;
724 
727  void removeInnerPrecName ();
728 
734  std::pair<Teuchos::ParameterList, bool> innerPrecParams () const;
735 
738  void removeInnerPrecParams ();
739 
741  static std::string defaultInnerPrecName ();
742 
747 
752 
754  Teuchos::RCP<row_matrix_type> ReorderedLocalizedMatrix_;
755 
757  Teuchos::RCP<row_matrix_type> innerMatrix_;
758 
760  bool IsInitialized_ = false;
762  bool IsComputed_ = false;
764  bool IsOverlapping_ = false;
766  int OverlapLevel_ = 0;
767 
774 
776  mutable Teuchos::RCP<const Teuchos::ParameterList> validParams_;
777 
783  Tpetra::CombineMode CombineMode_ = Tpetra::ZERO;
784  bool AvgOverlap_ = false;
786  bool UseReordering_ = false;
788  std::string ReorderingAlgorithm_ = "none";
790  bool FilterSingletons_ = false;
792  Teuchos::RCP<row_matrix_type> SingletonMatrix_;
794  int NumIterations_ = 1;
796  bool ZeroStartingSolution_ = true;
799 
801  int NumInitialize_ = 0;
803  int NumCompute_ = 0;
805  mutable int NumApply_ = 0;
807  double InitializeTime_ = 0.0;
809  double ComputeTime_ = 0.0;
811  mutable double ApplyTime_ = 0.0;
813  double InitializeFlops_ = 0.0;
815  double ComputeFlops_ = 0.0;
817  mutable double ApplyFlops_ = 0.0;
823  mutable std::unique_ptr<MV> overlapping_B_;
825  mutable std::unique_ptr<MV> overlapping_Y_;
827  mutable std::unique_ptr<MV> reduced_reordered_B_;
829  mutable std::unique_ptr<MV> reduced_reordered_Y_;
831  mutable std::unique_ptr<MV> reduced_B_;
833  mutable std::unique_ptr<MV> reduced_Y_;
835  mutable std::unique_ptr<MV> reordered_B_;
837  mutable std::unique_ptr<MV> reordered_Y_;
839  mutable std::unique_ptr<MV> num_overlap_copies_;
841  mutable std::unique_ptr<MV> R_;
843  mutable std::unique_ptr<MV> C_;
844 
851  mutable Teuchos::RCP<const import_type> DistributedImporter_;
852 }; // class AdditiveSchwarz
853 
854 }// end namespace
855 
856 #endif // IFPACK2_ADDITIVESCHWARZ_DECL_HPP
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
The Tpetra::RowMatrix specialization matching MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:299
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:60
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &plist)
Set the preconditioner&#39;s parameters.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:731
virtual void compute()
Computes all (coefficient) data necessary to apply the preconditioner.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1063
virtual ~AdditiveSchwarz()=default
Destructor.
virtual 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, putting the result in Y.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:292
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
The Tpetra::CrsMatrix specialization that is a subclass of MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:304
Declaration of interface for nested preconditioners.
virtual int getNumApply() const
Returns the number of calls to apply().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1151
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:36
virtual bool isInitialized() const
Returns true if the preconditioner has been successfully initialized, false otherwise.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1056
virtual double getInitializeTime() const
Returns the time spent in initialize().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1158
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1234
virtual Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const
The domain Map of this operator.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:230
Mix-in interface for nested preconditioners.
Definition: Ifpack2_Details_NestedPreconditioner.hpp:64
virtual double getComputeTime() const
Returns the time spent in compute().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1165
std::string description() const
Return a simple one-line description of this object.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1179
virtual std::ostream & print(std::ostream &os) const
Prints basic information on iostream. This function is used by operator&lt;&lt;.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1346
AdditiveSchwarz(const Teuchos::RCP< const row_matrix_type > &A)
Constructor that takes a matrix.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:215
virtual double getApplyTime() const
Returns the time spent in apply().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1172
typename MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:290
virtual void initialize()
Computes all (graph-related) data necessary to initialize the preconditioner.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:973
typename MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:287
virtual void setInnerPreconditioner(const Teuchos::RCP< Preconditioner< scalar_type, local_ordinal_type, global_ordinal_type, node_type > > &innerPrec)
Set the inner preconditioner.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1652
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Get a list of the preconditioner&#39;s default parameters.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:927
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:74
virtual Teuchos::RCP< const row_matrix_type > getMatrix() const
The input matrix.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:255
typename MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:284
Declaration of interface for preconditioners that can change their matrix after construction.
static const EVerbosityLevel verbLevel_default
typename Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:294
virtual Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const
The range Map of this operator.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:243
virtual void setParameters(const Teuchos::ParameterList &plist)
Set the preconditioner&#39;s parameters.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:718
virtual int getNumCompute() const
Returns the number of calls to compute().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1144
Additive Schwarz domain decomposition for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:250
virtual bool isComputed() const
Returns true if the preconditioner has been successfully computed, false otherwise.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1130
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1723
virtual int getOverlapLevel() const
Returns the level of overlap.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1356
typename MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:281
void setZeroStartingSolution(bool zeroStartingSolution)
Set this preconditioner&#39;s parameters.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:624
virtual int getNumInitialize() const
Returns the number of calls to initialize().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1137