43 #ifndef IFPACK2_BLOCKRELAXATION_DECL_HPP 
   44 #define IFPACK2_BLOCKRELAXATION_DECL_HPP 
   50 #include "Ifpack2_Partitioner.hpp" 
   52 #include "Ifpack2_ContainerFactory.hpp" 
   54 #include "Tpetra_BlockCrsMatrix.hpp" 
   55 #include <type_traits> 
   82 template<
class MatrixType, 
class ContainerType = Container<MatrixType> >
 
   85                                            typename MatrixType::local_ordinal_type,
 
   86                                            typename MatrixType::global_ordinal_type,
 
   87                                            typename MatrixType::node_type>,
 
   89                                                                        typename MatrixType::local_ordinal_type,
 
   90                                                                        typename MatrixType::global_ordinal_type,
 
   91                                                                        typename MatrixType::node_type> >
 
  112   typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> 
row_matrix_type;
 
  114   static_assert (std::is_same<MatrixType, row_matrix_type>::value,
 
  115                  "Ifpack2::BlockRelaxation: Please use MatrixType = Tpetra::RowMatrix.");
 
  118                  "Ifpack2::BlockRelaxation: Do NOT specify the (second) " 
  119                  "ContainerType template parameter explicitly.  The default " 
  120                  "value is fine.  Please instead specify the container type to " 
  121                  "use by setting the \"relaxation: container\" parameter.");
 
  124   typedef Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type> 
import_type;
 
  127   void computeImporter() 
const;
 
  133   typedef Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> vector_type;
 
  138   typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
 
  222     return(IsInitialized_);
 
  277   void apply(
const MV& X,
 
  289   bool hasTransposeApply() 
const;
 
  360   virtual void ApplyInverseJacobi (
const MV& X, MV& Y) 
const;
 
  362   virtual void ApplyInverseGS (
const MV& X, MV& Y) 
const;
 
  364   virtual void ApplyInverseSGS (
const MV& X, MV& Y) 
const;
 
  368   void ExtractSubmatricesStructure();
 
  395   std::string PartitionerType_;
 
  407   std::string containerType_;
 
  410   Details::RelaxationType PrecType_;
 
  416   bool ZeroStartingSolution_;
 
  420   bool hasBlockCrsMatrix_;
 
  447   mutable int NumApply_;
 
  450   double InitializeTime_;
 
  456   mutable double ApplyTime_;
 
  478 #endif // IFPACK2_BLOCKRELAXATION_DECL_HPP 
Mix-in interface for preconditioners that can change their matrix after construction. 
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
 
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const 
The communicator over which the input matrix is distributed. 
Definition: Ifpack2_BlockRelaxation_def.hpp:327
 
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:104
 
std::string description() const 
A one-line description of this object. 
Definition: Ifpack2_BlockRelaxation_def.hpp:963
 
double getApplyTime() const 
Returns the time spent in apply(). 
Definition: Ifpack2_BlockRelaxation_def.hpp:422
 
int getNumApply() const 
Returns the number of calls to apply(). 
Definition: Ifpack2_BlockRelaxation_def.hpp:398
 
Teuchos::RCP< const map_type > getRangeMap() const 
Returns the Tpetra::Map object associated with the range of this operator. 
Definition: Ifpack2_BlockRelaxation_def.hpp:364
 
int getNumCompute() const 
Returns the number of calls to compute(). 
Definition: Ifpack2_BlockRelaxation_def.hpp:390
 
void setParameters(const Teuchos::ParameterList ¶ms)
Sets all the parameters for the preconditioner. 
Definition: Ifpack2_BlockRelaxation_def.hpp:173
 
void compute()
compute the preconditioner for the specified matrix, diagonal perturbation thresholds and relaxation ...
Definition: Ifpack2_BlockRelaxation_def.hpp:662
 
Block relaxation preconditioners (or smoothers) for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse ma...
Definition: Ifpack2_BlockRelaxation_decl.hpp:83
 
double getInitializeTime() const 
Returns the time spent in initialize(). 
Definition: Ifpack2_BlockRelaxation_def.hpp:406
 
double getComputeTime() const 
Returns the time spent in compute(). 
Definition: Ifpack2_BlockRelaxation_def.hpp:414
 
bool isComputed() const 
Return true if compute() has been called. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:229
 
size_t getNodeSmootherComplexity() const 
Get a rough estimate of cost per iteration. 
Definition: Ifpack2_BlockRelaxation_def.hpp:429
 
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:109
 
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization corresponding to MatrixType. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:112
 
Tpetra::Import< local_ordinal_type, global_ordinal_type, node_type > import_type
Tpetra::Importer specialization for use with MatrixType and compatible MultiVectors. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:115
 
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned. 
Definition: Ifpack2_BlockRelaxation_def.hpp:57
 
virtual ~BlockRelaxation()
Destructor. 
Definition: Ifpack2_BlockRelaxation_def.hpp:118
 
Interface for all Ifpack2 preconditioners. 
Definition: Ifpack2_Preconditioner.hpp:107
 
Teuchos::RCP< const map_type > getDomainMap() const 
Returns the Tpetra::Map object associated with the domain of this operator. 
Definition: Ifpack2_BlockRelaxation_def.hpp:350
 
bool isInitialized() const 
Returns true if the preconditioner has been successfully initialized. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:221
 
int getNumInitialize() const 
Returns the number of calls to initialize(). 
Definition: Ifpack2_BlockRelaxation_def.hpp:383
 
Declaration of interface for preconditioners that can change their matrix after construction. 
 
void initialize()
Initialize. 
Definition: Ifpack2_BlockRelaxation_def.hpp:553
 
static const EVerbosityLevel verbLevel_default
 
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_BlockRelaxation_def.hpp:1020
 
MatrixType::node_type node_type
Node type of the input MatrixType. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:106
 
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:101
 
void applyMat(const MV &X, MV &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS) const 
Applies the matrix to a Tpetra::MultiVector. 
Definition: Ifpack2_BlockRelaxation_def.hpp:537
 
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:98
 
Interface for creating and solving a set of local linear problems. 
Definition: Ifpack2_Container_decl.hpp:112
 
Teuchos::RCP< const row_matrix_type > getMatrix() const 
The input matrix of this preconditioner's constructor. 
Definition: Ifpack2_BlockRelaxation_def.hpp:341
 
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const 
Return a list of all the parameters that this class accepts. 
Definition: Ifpack2_BlockRelaxation_def.hpp:124
 
void apply(const MV &X, MV &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const 
Applies the preconditioner to X, returns the result in Y. 
Definition: Ifpack2_BlockRelaxation_def.hpp:446
 
Teuchos::RCP< Ifpack2::Partitioner< Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > > > getPartitioner()
For diagnostic purposes. 
Definition: Ifpack2_BlockRelaxation_decl.hpp:349
 
BlockRelaxation(const Teuchos::RCP< const row_matrix_type > &Matrix)
Constructor. 
Definition: Ifpack2_BlockRelaxation_def.hpp:85