52 #ifndef IFPACK2_ADDITIVESCHWARZ_DECL_HPP 
   53 #define IFPACK2_ADDITIVESCHWARZ_DECL_HPP 
   55 #include "Ifpack2_ConfigDefs.hpp" 
   59 #include "Tpetra_Map.hpp" 
   60 #include "Tpetra_MultiVector.hpp" 
   61 #include "Tpetra_RowMatrix.hpp" 
   63 #include <type_traits> 
   67 template<
class MV, 
class OP, 
class NormType>
 
  276 template<
class MatrixType,
 
  277          class LocalInverseType =
 
  279                         typename MatrixType::local_ordinal_type,
 
  280                         typename MatrixType::global_ordinal_type,
 
  281                         typename MatrixType::node_type> >
 
  284                                   typename MatrixType::local_ordinal_type,
 
  285                                   typename MatrixType::global_ordinal_type,
 
  286                                   typename MatrixType::node_type>,
 
  288                                                               typename MatrixType::local_ordinal_type,
 
  289                                                               typename MatrixType::global_ordinal_type,
 
  290                                                               typename MatrixType::node_type> >,
 
  292                                                                 typename MatrixType::local_ordinal_type,
 
  293                                                                 typename MatrixType::global_ordinal_type,
 
  294                                                                 typename MatrixType::node_type> >
 
  297   static_assert(std::is_same<LocalInverseType,
 
  299                     typename MatrixType::local_ordinal_type,
 
  300                     typename MatrixType::global_ordinal_type,
 
  301                     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.");
 
  303   static_assert(std::is_same<MatrixType,
 
  304                   Tpetra::RowMatrix<
typename MatrixType::scalar_type,
 
  305                     typename MatrixType::local_ordinal_type,
 
  306                     typename MatrixType::global_ordinal_type,
 
  307                     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.");
 
  351                    const int overlapLevel);
 
  368   apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &X,
 
  369          Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> &Y,
 
  696   virtual std::ostream& 
print(std::ostream& os) 
const;
 
  737   void localApply(MV &OverlappingX, MV &OverlappingY) 
const;
 
  741   bool hasInnerPrecName () 
const;
 
  745   std::string innerPrecName () 
const;
 
  749   void removeInnerPrecName ();
 
  756   std::pair<Teuchos::ParameterList, bool> innerPrecParams () 
const;
 
  760   void removeInnerPrecParams ();
 
  763   static std::string defaultInnerPrecName ();
 
  782   bool IsInitialized_ = 
false;
 
  784   bool IsComputed_ = 
false;
 
  786   bool IsOverlapping_ = 
false;
 
  788   int OverlapLevel_ = 0;
 
  801   Tpetra::CombineMode CombineMode_ = Tpetra::ZERO;
 
  803   bool UseReordering_ = 
false;
 
  805   std::string ReorderingAlgorithm_ = 
"none";
 
  807   bool FilterSingletons_ = 
false;
 
  811   int NumIterations_ = 1;
 
  813   bool ZeroStartingSolution_ = 
true;
 
  818   int NumInitialize_ = 0;
 
  822   mutable int NumApply_ = 0;
 
  824   double InitializeTime_ = 0.0;
 
  826   double ComputeTime_ = 0.0;
 
  828   mutable double ApplyTime_ = 0.0;
 
  830   double InitializeFlops_ = 0.0;
 
  832   double ComputeFlops_ = 0.0;
 
  834   mutable double ApplyFlops_ = 0.0;
 
  840   mutable std::unique_ptr<MV> overlapping_B_;
 
  842   mutable std::unique_ptr<MV> overlapping_Y_;
 
  844   mutable std::unique_ptr<MV> R_;
 
  846   mutable std::unique_ptr<MV> C_;
 
  859 #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:331
 
Mix-in interface for preconditioners that can change their matrix after construction. 
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
 
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &plist)
Set the preconditioner's parameters. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:711
 
virtual void compute()
Computes all (coefficient) data necessary to apply the preconditioner. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1009
 
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:294
 
Declaration of interface for nested preconditioners. 
 
virtual int getNumApply() const 
Returns the number of calls to apply(). 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1079
 
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:68
 
virtual bool isInitialized() const 
Returns true if the preconditioner has been successfully initialized, false otherwise. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1002
 
virtual double getInitializeTime() const 
Returns the time spent in initialize(). 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1086
 
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:1162
 
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:260
 
Mix-in interface for nested preconditioners. 
Definition: Ifpack2_Details_NestedPreconditioner.hpp:97
 
virtual double getComputeTime() const 
Returns the time spent in compute(). 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1093
 
std::string description() const 
Return a simple one-line description of this object. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1107
 
virtual std::ostream & print(std::ostream &os) const 
Prints basic information on iostream. This function is used by operator<<. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1269
 
AdditiveSchwarz(const Teuchos::RCP< const row_matrix_type > &A)
Constructor that takes a matrix. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:245
 
virtual double getApplyTime() const 
Returns the time spent in apply(). 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1100
 
typename MatrixType::node_type node_type
The Node type used by the input MatrixType. 
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:322
 
virtual void initialize()
Computes all (graph-related) data necessary to initialize the preconditioner. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:924
 
typename MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType. 
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:319
 
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:1479
 
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const 
Get a list of the preconditioner's default parameters. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:878
 
Interface for all Ifpack2 preconditioners. 
Definition: Ifpack2_Preconditioner.hpp:107
 
virtual Teuchos::RCP< const row_matrix_type > getMatrix() const 
The input matrix. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:285
 
typename MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType. 
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:316
 
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:326
 
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:273
 
virtual void setParameters(const Teuchos::ParameterList &plist)
Set the preconditioner's parameters. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:698
 
virtual int getNumCompute() const 
Returns the number of calls to compute(). 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1072
 
Additive Schwarz domain decomposition for Tpetra sparse matrices. 
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:282
 
virtual bool isComputed() const 
Returns true if the preconditioner has been successfully computed, false otherwise. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1058
 
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1547
 
virtual int getOverlapLevel() const 
Returns the level of overlap. 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1279
 
typename MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType. 
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:313
 
virtual int getNumInitialize() const 
Returns the number of calls to initialize(). 
Definition: Ifpack2_AdditiveSchwarz_def.hpp:1065