10#ifndef IFPACK2_BLOCKRELAXATION_DECL_HPP
11#define IFPACK2_BLOCKRELAXATION_DECL_HPP
17#include "Ifpack2_Partitioner.hpp"
19#include "Ifpack2_ContainerFactory.hpp"
20#include "Tpetra_BlockCrsMatrix.hpp"
48template <
class MatrixType,
class ContainerType = Container<MatrixType>>
50 typename MatrixType::local_ordinal_type,
51 typename MatrixType::global_ordinal_type,
52 typename MatrixType::node_type>,
54 typename MatrixType::local_ordinal_type,
55 typename MatrixType::global_ordinal_type,
56 typename MatrixType::node_type>> {
73 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType
magnitude_type;
76 typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type>
row_matrix_type;
78 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
79 "Ifpack2::BlockRelaxation: Please use MatrixType = Tpetra::RowMatrix.");
81 static_assert(std::is_same<ContainerType, Container<row_matrix_type>>::value,
82 "Ifpack2::BlockRelaxation: Do NOT specify the (second) "
83 "ContainerType template parameter explicitly. The default "
84 "value is fine. Please instead specify the container type to "
85 "use by setting the \"relaxation: container\" parameter.");
88 typedef Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type>
import_type;
95 void setParametersImpl(Teuchos::ParameterList& params);
97 void computeImporter()
const;
104 typedef Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> vector_type;
105 typedef Teuchos::ScalarTraits<scalar_type> STS;
106 typedef Teuchos::ScalarTraits<magnitude_type> STM;
112 block_crs_matrix_type;
113 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
150 explicit BlockRelaxation(
const Teuchos::RCP<const row_matrix_type>& Matrix);
189 bool supportsZeroStartingSolution() {
return true; }
194 Teuchos::RCP<const Teuchos::ParameterList>
202 return (IsInitialized_);
210 return (IsComputed_);
240 setMatrix(
const Teuchos::RCP<const row_matrix_type>& A);
257 void apply(
const MV& X,
259 Teuchos::ETransp mode = Teuchos::NO_TRANS,
260 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
261 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero())
const;
269 bool hasTransposeApply()
const;
280 Teuchos::ETransp mode = Teuchos::NO_TRANS)
const;
287 Teuchos::RCP<const Teuchos::Comm<int>>
getComm()
const;
290 Teuchos::RCP<const row_matrix_type>
getMatrix()
const;
322 describe(Teuchos::FancyOStream& out,
323 const Teuchos::EVerbosityLevel verbLevel =
324 Teuchos::Describable::verbLevel_default)
const;
329 Teuchos::RCP<Ifpack2::Partitioner<Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>>>
getPartitioner() {
return Partitioner_; }
339 virtual void ApplyInverseJacobi(
const MV& X, MV& Y)
const;
341 virtual void ApplyInverseGS(
const MV& X, MV& Y)
const;
343 virtual void ApplyInverseSGS(
const MV& X, MV& Y)
const;
347 void ExtractSubmatricesStructure();
354 Teuchos::RCP<const row_matrix_type> A_;
357 mutable Teuchos::RCP<Container<row_matrix_type>> Container_;
364 mutable Teuchos::RCP<vector_type> DiagRCP_;
367 Teuchos::RCP<Ifpack2::Partitioner<Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type>>> Partitioner_;
371 std::string PartitionerType_;
374 Teuchos::ParameterList List_;
383 std::string containerType_;
386 Details::RelaxationType PrecType_;
392 bool ZeroStartingSolution_;
396 bool hasBlockCrsMatrix_;
409 std::string schwarzCombineMode_;
433 mutable int NumApply_;
436 double InitializeTime_;
442 mutable double ApplyTime_;
455 Teuchos::RCP<vector_type> W_;
457 mutable Teuchos::RCP<const Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type>> Importer_;
Declaration of interface for preconditioners that can change their matrix after construction.
Block relaxation preconditioners (or smoothers) for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse ma...
Definition Ifpack2_BlockRelaxation_decl.hpp:56
Teuchos::RCP< const row_matrix_type > getMatrix() const
The input matrix of this preconditioner's constructor.
Definition Ifpack2_BlockRelaxation_def.hpp:340
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:76
Teuchos::RCP< const map_type > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition Ifpack2_BlockRelaxation_def.hpp:349
void setZeroStartingSolution(bool zeroStartingSolution)
Set this preconditioner's parameters.
Definition Ifpack2_BlockRelaxation_decl.hpp:191
int getNumApply() const
Returns the number of calls to apply().
Definition Ifpack2_BlockRelaxation_def.hpp:390
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition Ifpack2_BlockRelaxation_def.hpp:96
void compute()
compute the preconditioner for the specified matrix, diagonal perturbation thresholds and relaxation ...
Definition Ifpack2_BlockRelaxation_def.hpp:742
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition Ifpack2_BlockRelaxation_def.hpp:416
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:1068
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:88
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:65
virtual ~BlockRelaxation()
Destructor.
Definition Ifpack2_BlockRelaxation_def.hpp:91
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition Ifpack2_BlockRelaxation_decl.hpp:201
int getNumCompute() const
Returns the number of calls to compute().
Definition Ifpack2_BlockRelaxation_def.hpp:384
double getApplyTime() const
Returns the time spent in apply().
Definition Ifpack2_BlockRelaxation_def.hpp:411
int getNumInitialize() const
Returns the number of calls to initialize().
Definition Ifpack2_BlockRelaxation_def.hpp:378
MatrixType::node_type node_type
Node type of the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:70
void initialize()
Initialize.
Definition Ifpack2_BlockRelaxation_def.hpp:555
Teuchos::RCP< const map_type > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition Ifpack2_BlockRelaxation_def.hpp:362
bool isComputed() const
Return true if compute() has been called.
Definition Ifpack2_BlockRelaxation_decl.hpp:209
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:541
double getComputeTime() const
Returns the time spent in compute().
Definition Ifpack2_BlockRelaxation_def.hpp:404
void setParameters(const Teuchos::ParameterList ¶ms)
Sets all the parameters for the preconditioner.
Definition Ifpack2_BlockRelaxation_def.hpp:156
double getInitializeTime() const
Returns the time spent in initialize().
Definition Ifpack2_BlockRelaxation_def.hpp:397
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:433
std::string description() const
A one-line description of this object.
Definition Ifpack2_BlockRelaxation_def.hpp:1007
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:68
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the input matrix is distributed.
Definition Ifpack2_BlockRelaxation_def.hpp:327
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition Ifpack2_BlockRelaxation_decl.hpp:73
Teuchos::RCP< Ifpack2::Partitioner< Tpetra::RowGraph< local_ordinal_type, global_ordinal_type, node_type > > > getPartitioner()
For diagnostic purposes.
Definition Ifpack2_BlockRelaxation_decl.hpp:329
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_BlockRelaxation_def.hpp:31
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition Ifpack2_BlockRelaxation_decl.hpp:62
Mix-in interface for preconditioners that can change their matrix after construction.
Definition Ifpack2_Details_CanChangeMatrix.hpp:60
Interface for all Ifpack2 preconditioners.
Definition Ifpack2_Preconditioner.hpp:74
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40