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 namespace Ifpack2 {
41 
243 template <class MatrixType,
244  class LocalInverseType =
245  Preconditioner<typename MatrixType::scalar_type,
246  typename MatrixType::local_ordinal_type,
247  typename MatrixType::global_ordinal_type,
248  typename MatrixType::node_type>>
249 class AdditiveSchwarz : virtual public Preconditioner<typename MatrixType::scalar_type,
250  typename MatrixType::local_ordinal_type,
251  typename MatrixType::global_ordinal_type,
252  typename MatrixType::node_type>,
253  virtual public Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
254  typename MatrixType::local_ordinal_type,
255  typename MatrixType::global_ordinal_type,
256  typename MatrixType::node_type>>,
257  virtual public Details::NestedPreconditioner<Preconditioner<typename MatrixType::scalar_type,
258  typename MatrixType::local_ordinal_type,
259  typename MatrixType::global_ordinal_type,
260  typename MatrixType::node_type>> {
261  public:
262  static_assert(std::is_same<LocalInverseType,
263  Preconditioner<typename MatrixType::scalar_type,
264  typename MatrixType::local_ordinal_type,
265  typename MatrixType::global_ordinal_type,
266  typename MatrixType::node_type>>::value,
267  "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.");
268 
269  static_assert(std::is_same<MatrixType,
270  Tpetra::RowMatrix<typename MatrixType::scalar_type,
271  typename MatrixType::local_ordinal_type,
272  typename MatrixType::global_ordinal_type,
273  typename MatrixType::node_type>>::value,
274  "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.");
275 
277 
278 
280  using scalar_type = typename MatrixType::scalar_type;
281 
283  using local_ordinal_type = typename MatrixType::local_ordinal_type;
284 
286  using global_ordinal_type = typename MatrixType::global_ordinal_type;
287 
289  using node_type = typename MatrixType::node_type;
290 
292  using magnitude_type =
294 
296  using row_matrix_type =
297  Tpetra::RowMatrix<scalar_type, local_ordinal_type,
299 
301  using crs_matrix_type =
302  Tpetra::CrsMatrix<scalar_type, local_ordinal_type,
304 
306  using coord_type =
307  Tpetra::MultiVector<magnitude_type, local_ordinal_type,
309 
311  // \name Constructors and destructor
313 
318 
324  const Teuchos::RCP<const coord_type>& coordinates);
325 
336  const int overlapLevel);
337 
339  virtual ~AdditiveSchwarz() = default;
340 
342 
344 
347 
350 
352  virtual void
353  apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
354  Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
358 
360 
362 
397  virtual void
401  node_type>>& innerPrec);
402 
404 
406 
429  virtual void
432 
435 
439  void
440  setCoord(const Teuchos::RCP<const coord_type>& Coordinates);
441 
444 
614  virtual void setParameters(const Teuchos::ParameterList& plist);
615 
639  void
641 
642  bool supportsZeroStartingSolution() { return true; }
643 
644  void setZeroStartingSolution(bool zeroStartingSolution) { ZeroStartingSolution_ = zeroStartingSolution; };
645 
651 
653  virtual void initialize();
654 
656  virtual bool isInitialized() const;
657 
659  virtual void compute();
660 
662  virtual bool isComputed() const;
663 
665  virtual int getNumInitialize() const;
666 
668  virtual int getNumCompute() const;
669 
671  virtual int getNumApply() const;
672 
674  virtual double getInitializeTime() const;
675 
677  virtual double getComputeTime() const;
678 
680  virtual double getApplyTime() const;
681 
683 
684 
686  std::string description() const;
687 
690 
692 
694  virtual std::ostream& print(std::ostream& os) const;
695 
697  virtual int getOverlapLevel() const;
698 
699  private:
701  typedef Tpetra::Map<local_ordinal_type,
703  node_type>
704  map_type;
706  typedef Tpetra::Import<local_ordinal_type,
708  node_type>
709  import_type;
711  typedef Tpetra::MultiVector<scalar_type,
714  node_type>
715  MV;
717  typedef Tpetra::Operator<scalar_type,
720  node_type>
721  OP;
723  typedef Preconditioner<scalar_type,
726  node_type>
727  prec_type;
728 
731 
733  typedef Kokkos::DualView<local_ordinal_type*, typename crs_matrix_type::device_type> perm_dualview_type;
734 
736  AdditiveSchwarz(const AdditiveSchwarz& RHS);
737 
739  void setup();
740 
742  void localApply(MV& OverlappingX, MV& OverlappingY) const;
743 
746  bool hasInnerPrecName() const;
747 
749  std::string innerPrecName() const;
750 
753  void removeInnerPrecName();
754 
760  std::pair<Teuchos::ParameterList, bool> innerPrecParams() const;
761 
764  void removeInnerPrecParams();
765 
767  static std::string defaultInnerPrecName();
768 
773 
778 
780  Teuchos::RCP<row_matrix_type> ReorderedLocalizedMatrix_;
781 
783  Teuchos::RCP<row_matrix_type> innerMatrix_;
784 
787  Teuchos::RCP<const coord_type> Coordinates_ = Teuchos::null;
788 
790  bool IsInitialized_ = false;
792  bool IsComputed_ = false;
794  bool IsOverlapping_ = false;
796  int OverlapLevel_ = 0;
797 
804 
806  mutable Teuchos::RCP<const Teuchos::ParameterList> validParams_;
807 
813  Tpetra::CombineMode CombineMode_ = Tpetra::ZERO;
814  bool AvgOverlap_ = false;
816  bool UseReordering_ = false;
818  std::string ReorderingAlgorithm_ = "none";
820  bool FilterSingletons_ = false;
822  Teuchos::RCP<row_matrix_type> SingletonMatrix_;
824  int NumIterations_ = 1;
826  bool ZeroStartingSolution_ = true;
829 
831  int NumInitialize_ = 0;
833  int NumCompute_ = 0;
835  mutable int NumApply_ = 0;
837  double InitializeTime_ = 0.0;
839  double ComputeTime_ = 0.0;
841  mutable double ApplyTime_ = 0.0;
843  double InitializeFlops_ = 0.0;
845  double ComputeFlops_ = 0.0;
847  mutable double ApplyFlops_ = 0.0;
853  mutable std::unique_ptr<MV> overlapping_B_;
855  mutable std::unique_ptr<MV> overlapping_Y_;
857  mutable std::unique_ptr<MV> reduced_reordered_B_;
859  mutable std::unique_ptr<MV> reduced_reordered_Y_;
861  mutable std::unique_ptr<MV> reduced_B_;
863  mutable std::unique_ptr<MV> reduced_Y_;
865  mutable std::unique_ptr<MV> reordered_B_;
867  mutable std::unique_ptr<MV> reordered_Y_;
869  mutable std::unique_ptr<MV> num_overlap_copies_;
871  mutable std::unique_ptr<MV> R_;
873  mutable std::unique_ptr<MV> C_;
875  perm_dualview_type perm_coors;
876 
883  mutable Teuchos::RCP<const import_type> DistributedImporter_;
884 }; // class AdditiveSchwarz
885 
886 } // namespace Ifpack2
887 
888 #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:298
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:690
virtual void compute()
Computes all (coefficient) data necessary to apply the preconditioner.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1051
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:273
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:303
Declaration of interface for nested preconditioners.
virtual int getNumApply() const
Returns the number of calls to apply().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1132
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:1046
virtual double getInitializeTime() const
Returns the time spent in initialize().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1137
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:1205
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:214
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:1142
std::string description() const
Return a simple one-line description of this object.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1152
Tpetra::MultiVector< magnitude_type, local_ordinal_type, global_ordinal_type, node_type > coord_type
The Tpetra::MultiVector specialization used for containing coordinates.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:308
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:1313
AdditiveSchwarz(const Teuchos::RCP< const row_matrix_type > &A)
Constructor that takes a matrix.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:194
virtual double getApplyTime() const
Returns the time spent in apply().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1147
typename MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:289
virtual void initialize()
Computes all (graph-related) data necessary to initialize the preconditioner.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:925
typename MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:286
Teuchos::RCP< const coord_type > getCoord() const
Get the coordinates associated with the input matrix&#39;s rows.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:242
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Get a list of the preconditioner&#39;s default parameters.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:881
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:237
typename MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:283
Declaration of interface for preconditioners that can change their matrix after construction.
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:1618
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:293
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:226
virtual void setParameters(const Teuchos::ParameterList &plist)
Set the preconditioner&#39;s parameters.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:680
virtual int getNumCompute() const
Returns the number of calls to compute().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1127
Additive Schwarz domain decomposition for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:249
virtual bool isComputed() const
Returns true if the preconditioner has been successfully computed, false otherwise.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1117
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1690
virtual int getOverlapLevel() const
Returns the level of overlap.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1321
void setCoord(const Teuchos::RCP< const coord_type > &Coordinates)
Set the matrix rows&#39; coordinates.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1714
typename MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:280
void setZeroStartingSolution(bool zeroStartingSolution)
Set this preconditioner&#39;s parameters.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:644
virtual int getNumInitialize() const
Returns the number of calls to initialize().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1122