47 #ifndef IFPACK2_HIPTMAIR_DECL_HPP 
   48 #define IFPACK2_HIPTMAIR_DECL_HPP 
   51 #include <type_traits> 
   69   template<
class MatrixType>
 
   72                                            typename MatrixType::local_ordinal_type,
 
   73                                            typename MatrixType::global_ordinal_type,
 
   74                                            typename MatrixType::node_type>
 
  101     static_assert(std::is_same<MatrixType, row_matrix_type>::value, 
"Ifpack2::Hiptmair: The template parameter MatrixType must be a Tpetra::RowMatrix specialization.  Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore.  The constructor can take either a RowMatrix or a CrsMatrix just fine.");
 
  141       return IsInitialized_;
 
  158     apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
 
  159            Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
 
  165     applyHiptmairSmoother(
const Tpetra::MultiVector<
typename MatrixType::scalar_type,
 
  166                           typename MatrixType::local_ordinal_type,
 
  167                           typename MatrixType::global_ordinal_type,
 
  168                           typename MatrixType::node_type>& X,
 
  169                           Tpetra::MultiVector<
typename MatrixType::scalar_type,
 
  170                           typename MatrixType::local_ordinal_type,
 
  171                           typename MatrixType::global_ordinal_type,
 
  172                           typename MatrixType::node_type>& Y) 
const;
 
  239     std::string precType1_, precType2_, preOrPost_;
 
  242     bool ZeroStartingSolution_;
 
  256     mutable int NumApply_;
 
  258     double InitializeTime_;
 
  262     mutable double ApplyTime_;
 
  269 #endif // IFPACK2_HIPTMAIR_DECL_HPP 
void initialize()
Do any initialization that depends on the input matrix's structure. 
Definition: Ifpack2_Hiptmair_def.hpp:210
 
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType. 
Definition: Ifpack2_Hiptmair_decl.hpp:84
 
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType. 
Definition: Ifpack2_Hiptmair_decl.hpp:81
 
double getComputeTime() const 
Returns the time spent in Compute(). 
Definition: Ifpack2_Hiptmair_def.hpp:198
 
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const 
Returns the operator's communicator. 
Definition: Ifpack2_Hiptmair_def.hpp:125
 
bool hasTransposeApply() const 
Whether this object's apply() method can apply the transpose (or conjugate transpose, if applicable). 
Definition: Ifpack2_Hiptmair_def.hpp:166
 
bool isInitialized() const 
Return true if initialize() completed successfully, else false. 
Definition: Ifpack2_Hiptmair_decl.hpp:140
 
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_Hiptmair_def.hpp:428
 
double getApplyTime() const 
Returns the time spent in apply(). 
Definition: Ifpack2_Hiptmair_def.hpp:204
 
void setParameters(const Teuchos::ParameterList ¶ms)
Set the preconditioner's parameters. 
Definition: Ifpack2_Hiptmair_def.hpp:84
 
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Type of the Tpetra::RowMatrix specialization that this class uses. 
Definition: Ifpack2_Hiptmair_decl.hpp:99
 
Teuchos::RCP< const Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > > getMatrix() const 
Returns a reference to the matrix to be preconditioned. 
Definition: Ifpack2_Hiptmair_def.hpp:136
 
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const 
Tpetra::Map representing the range of this operator. 
Definition: Ifpack2_Hiptmair_def.hpp:155
 
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry. 
Definition: Ifpack2_Hiptmair_decl.hpp:93
 
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType. 
Definition: Ifpack2_Hiptmair_decl.hpp:87
 
int getNumApply() const 
Returns the number of calls to apply(). 
Definition: Ifpack2_Hiptmair_def.hpp:186
 
Wrapper for Hiptmair smoothers. 
Definition: Ifpack2_Hiptmair_decl.hpp:70
 
Hiptmair(const Teuchos::RCP< const row_matrix_type > &A, const Teuchos::RCP< const row_matrix_type > &PtAP, const Teuchos::RCP< const row_matrix_type > &P)
Constructor that takes 3 Tpetra matrices. 
Definition: Ifpack2_Hiptmair_def.hpp:57
 
int getNumCompute() const 
Returns the number of calls to Compute(). 
Definition: Ifpack2_Hiptmair_def.hpp:180
 
int getNumInitialize() const 
Returns the number of calls to Initialize(). 
Definition: Ifpack2_Hiptmair_def.hpp:174
 
MatrixType::node_type node_type
The Node type used by the input MatrixType. 
Definition: Ifpack2_Hiptmair_decl.hpp:90
 
Interface for all Ifpack2 preconditioners. 
Definition: Ifpack2_Preconditioner.hpp:107
 
bool isComputed() const 
Return true if compute() completed successfully, else false. 
Definition: Ifpack2_Hiptmair_decl.hpp:148
 
void compute()
Do any initialization that depends on the input matrix's values. 
Definition: Ifpack2_Hiptmair_def.hpp:247
 
static const EVerbosityLevel verbLevel_default
 
virtual ~Hiptmair()
Destructor. 
Definition: Ifpack2_Hiptmair_def.hpp:81
 
std::string description() const 
Return a simple one-line description of this object. 
Definition: Ifpack2_Hiptmair_def.hpp:398
 
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const 
Tpetra::Map representing the domain of this operator. 
Definition: Ifpack2_Hiptmair_def.hpp:143
 
Ifpack2::Preconditioner< scalar_type, local_ordinal_type, global_ordinal_type, node_type > prec_type
Type of the Ifpack2::Preconditioner specialization from which this class inherits. 
Definition: Ifpack2_Hiptmair_decl.hpp:101
 
double getInitializeTime() const 
Returns the time spent in Initialize(). 
Definition: Ifpack2_Hiptmair_def.hpp:192
 
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_Hiptmair_def.hpp:273