19 #ifndef IFPACK2_ADDITIVESCHWARZ_DECL_HPP
20 #define IFPACK2_ADDITIVESCHWARZ_DECL_HPP
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"
31 #include <type_traits>
35 template <
class MV,
class OP,
class NormType>
243 template <
class MatrixType,
244 class LocalInverseType =
246 typename MatrixType::local_ordinal_type,
247 typename MatrixType::global_ordinal_type,
248 typename MatrixType::node_type>>
250 typename MatrixType::local_ordinal_type,
251 typename MatrixType::global_ordinal_type,
252 typename MatrixType::node_type>,
254 typename MatrixType::local_ordinal_type,
255 typename MatrixType::global_ordinal_type,
256 typename MatrixType::node_type>>,
258 typename MatrixType::local_ordinal_type,
259 typename MatrixType::global_ordinal_type,
260 typename MatrixType::node_type>> {
262 static_assert(std::is_same<LocalInverseType,
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.");
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.");
324 const int overlapLevel);
341 apply(
const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
342 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
621 bool supportsZeroStartingSolution() {
return true; }
673 virtual std::ostream&
print(std::ostream& os)
const;
718 void localApply(MV& OverlappingX, MV& OverlappingY)
const;
722 bool hasInnerPrecName()
const;
725 std::string innerPrecName()
const;
729 void removeInnerPrecName();
736 std::pair<Teuchos::ParameterList, bool> innerPrecParams()
const;
740 void removeInnerPrecParams();
743 static std::string defaultInnerPrecName();
762 bool IsInitialized_ =
false;
764 bool IsComputed_ =
false;
766 bool IsOverlapping_ =
false;
768 int OverlapLevel_ = 0;
785 Tpetra::CombineMode CombineMode_ = Tpetra::ZERO;
786 bool AvgOverlap_ =
false;
788 bool UseReordering_ =
false;
790 std::string ReorderingAlgorithm_ =
"none";
792 bool FilterSingletons_ =
false;
796 int NumIterations_ = 1;
798 bool ZeroStartingSolution_ =
true;
803 int NumInitialize_ = 0;
807 mutable int NumApply_ = 0;
809 double InitializeTime_ = 0.0;
811 double ComputeTime_ = 0.0;
813 mutable double ApplyTime_ = 0.0;
815 double InitializeFlops_ = 0.0;
817 double ComputeFlops_ = 0.0;
819 mutable double ApplyFlops_ = 0.0;
825 mutable std::unique_ptr<MV> overlapping_B_;
827 mutable std::unique_ptr<MV> overlapping_Y_;
829 mutable std::unique_ptr<MV> reduced_reordered_B_;
831 mutable std::unique_ptr<MV> reduced_reordered_Y_;
833 mutable std::unique_ptr<MV> reduced_B_;
835 mutable std::unique_ptr<MV> reduced_Y_;
837 mutable std::unique_ptr<MV> reordered_B_;
839 mutable std::unique_ptr<MV> reordered_Y_;
841 mutable std::unique_ptr<MV> num_overlap_copies_;
843 mutable std::unique_ptr<MV> R_;
845 mutable std::unique_ptr<MV> C_;
858 #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's parameters.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:678
virtual void compute()
Computes all (coefficient) data necessary to apply the preconditioner.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1001
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:261
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:1082
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:996
virtual double getInitializeTime() const
Returns the time spent in initialize().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1087
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:1155
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:207
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:1092
std::string description() const
Return a simple one-line description of this object.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1102
virtual std::ostream & print(std::ostream &os) const
Prints basic information on iostream. This function is used by operator<<.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1263
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:1097
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:913
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 Teuchos::ParameterList > getValidParameters() const
Get a list of the preconditioner's default parameters.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:869
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:230
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:1557
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:219
virtual void setParameters(const Teuchos::ParameterList &plist)
Set the preconditioner's parameters.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:668
virtual int getNumCompute() const
Returns the number of calls to compute().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1077
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:1067
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1629
virtual int getOverlapLevel() const
Returns the level of overlap.
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1271
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's parameters.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:623
virtual int getNumInitialize() const
Returns the number of calls to initialize().
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1072