Ifpack2 Templated Preconditioning Package
Version 1.0
|
Wrapper for Hiptmair smoothers. More...
#include <Ifpack2_Hiptmair_decl.hpp>
Public Types | |
typedef MatrixType::scalar_type | scalar_type |
The type of the entries of the input MatrixType. More... | |
typedef MatrixType::local_ordinal_type | local_ordinal_type |
The type of local indices in the input MatrixType. More... | |
typedef MatrixType::global_ordinal_type | global_ordinal_type |
The type of global indices in the input MatrixType. More... | |
typedef MatrixType::node_type | node_type |
The Node type used by the input MatrixType. More... | |
typedef Teuchos::ScalarTraits < scalar_type >::magnitudeType | magnitude_type |
The type of the magnitude (absolute value) of a matrix entry. More... | |
typedef 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. More... | |
typedef 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. More... | |
Public Types inherited from Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > | |
typedef Teuchos::ScalarTraits < MatrixType::scalar_type > ::magnitudeType | magnitude_type |
The type of the magnitude (absolute value) of a matrix entry. More... | |
Public Member Functions | |
Hiptmair (const Teuchos::RCP< const row_matrix_type > &A) | |
Constructor that takes 1 Tpetra matrix (assumes we'll get the rest off the parameter list) More... | |
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, const Teuchos::RCP< const row_matrix_type > &Pt=Teuchos::null) | |
Constructor that takes 3 Tpetra matrices. More... | |
virtual | ~Hiptmair () |
Destructor. More... | |
Methods for setting parameters and initialization | |
void | setParameters (const Teuchos::ParameterList ¶ms) |
Set the preconditioner's parameters. More... | |
bool | supportsZeroStartingSolution () |
void | setZeroStartingSolution (bool zeroStartingSolution) |
Set this preconditioner's parameters. More... | |
void | initialize () |
Do any initialization that depends on the input matrix's structure. More... | |
bool | isInitialized () const |
Return true if initialize() completed successfully, else false . More... | |
void | compute () |
Do any initialization that depends on the input matrix's values. More... | |
bool | isComputed () const |
Return true if compute() completed successfully, else false . More... | |
Implementation of Tpetra::Operator | |
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. More... | |
void | applyHiptmairSmoother (const Tpetra::MultiVector< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > &X, Tpetra::MultiVector< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > &Y) const |
void | updateCachedMultiVectors (const Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > &map1, const Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type >> &map2, size_t numVecs) const |
A service routine for updating the cached MultiVectors. More... | |
Teuchos::RCP< const Tpetra::Map < local_ordinal_type, global_ordinal_type, node_type > > | getDomainMap () const |
Tpetra::Map representing the domain of this operator. More... | |
Teuchos::RCP< const Tpetra::Map < local_ordinal_type, global_ordinal_type, node_type > > | getRangeMap () const |
Tpetra::Map representing the range of this operator. More... | |
bool | hasTransposeApply () const |
Whether this object's apply() method can apply the transpose (or conjugate transpose, if applicable). More... | |
Mathematical functions. | |
Teuchos::RCP< const Teuchos::Comm< int > > | getComm () const |
Returns the operator's communicator. More... | |
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. More... | |
int | getNumInitialize () const |
Returns the number of calls to Initialize(). More... | |
int | getNumCompute () const |
Returns the number of calls to Compute(). More... | |
int | getNumApply () const |
Returns the number of calls to apply(). More... | |
double | getInitializeTime () const |
Returns the time spent in Initialize(). More... | |
double | getComputeTime () const |
Returns the time spent in Compute(). More... | |
double | getApplyTime () const |
Returns the time spent in apply(). More... | |
Teuchos::RCP < Ifpack2::Preconditioner < typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > | getPrec1 () |
Returns prec 1. More... | |
Teuchos::RCP < Ifpack2::Preconditioner < typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > | getPrec2 () |
Returns prec 2. More... | |
Overridden from Teuchos::Describable | |
std::string | description () const |
Return a simple one-line description of this object. More... | |
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. More... | |
Public Member Functions inherited from Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > | |
virtual | ~Preconditioner () |
Destructor. More... | |
Wrapper for Hiptmair smoothers.
A | specialization of Tpetra::RowMatrix. |
Ifpack2::Hiptmair does smoothing on two spaces; a primary space and an auxiliary space. This situation arises when preconditioning Maxwell's equations discretized by edge elements.
For a list of all run-time parameters that this class accepts, see the documentation of setParameters().
typedef MatrixType::scalar_type Ifpack2::Hiptmair< MatrixType >::scalar_type |
The type of the entries of the input MatrixType.
typedef MatrixType::local_ordinal_type Ifpack2::Hiptmair< MatrixType >::local_ordinal_type |
The type of local indices in the input MatrixType.
typedef MatrixType::global_ordinal_type Ifpack2::Hiptmair< MatrixType >::global_ordinal_type |
The type of global indices in the input MatrixType.
typedef MatrixType::node_type Ifpack2::Hiptmair< MatrixType >::node_type |
The Node type used by the input MatrixType.
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Ifpack2::Hiptmair< MatrixType >::magnitude_type |
The type of the magnitude (absolute value) of a matrix entry.
typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::Hiptmair< MatrixType >::row_matrix_type |
Type of the Tpetra::RowMatrix specialization that this class uses.
typedef Ifpack2::Preconditioner<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::Hiptmair< MatrixType >::prec_type |
Type of the Ifpack2::Preconditioner specialization from which this class inherits.
|
explicit |
Constructor that takes 1 Tpetra matrix (assumes we'll get the rest off the parameter list)
|
explicit |
Constructor that takes 3 Tpetra matrices.
|
virtual |
Destructor.
|
virtual |
Set the preconditioner's parameters.
This preconditioner accepts the following parameters:
std::string
)std::string
)Teuchos::ParameterList
)Teuchos::ParameterList
)std::string
)bool
)
|
inlinevirtual |
Set this preconditioner's parameters.
Reimplemented from Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >.
|
virtual |
Do any initialization that depends on the input matrix's structure.
|
inlinevirtual |
Return true
if initialize() completed successfully, else false
.
|
virtual |
Do any initialization that depends on the input matrix's values.
|
inlinevirtual |
Return true
if compute() completed successfully, else false
.
|
virtual |
Apply the preconditioner to X, putting the result in Y.
void Ifpack2::Hiptmair< MatrixType >::updateCachedMultiVectors | ( | const Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > & | map1, |
const Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type >> & | map2, | ||
size_t | numVecs | ||
) | const |
A service routine for updating the cached MultiVectors.
|
virtual |
Tpetra::Map representing the domain of this operator.
|
virtual |
Tpetra::Map representing the range of this operator.
bool Ifpack2::Hiptmair< MatrixType >::hasTransposeApply | ( | ) | const |
Whether this object's apply() method can apply the transpose (or conjugate transpose, if applicable).
Teuchos::RCP< const Teuchos::Comm< int > > Ifpack2::Hiptmair< MatrixType >::getComm | ( | ) | const |
Returns the operator's communicator.
|
virtual |
Returns a reference to the matrix to be preconditioned.
|
virtual |
Returns the number of calls to Initialize().
|
virtual |
Returns the number of calls to Compute().
|
virtual |
Returns the number of calls to apply().
|
virtual |
Returns the time spent in Initialize().
|
virtual |
Returns the time spent in Compute().
|
virtual |
Returns the time spent in apply().
Teuchos::RCP< Ifpack2::Preconditioner< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::Hiptmair< MatrixType >::getPrec1 | ( | ) |
Returns prec 1.
Teuchos::RCP< Ifpack2::Preconditioner< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::Hiptmair< MatrixType >::getPrec2 | ( | ) |
Returns prec 2.
std::string Ifpack2::Hiptmair< MatrixType >::description | ( | ) | const |
Return a simple one-line description of this object.
void Ifpack2::Hiptmair< MatrixType >::describe | ( | Teuchos::FancyOStream & | out, |
const Teuchos::EVerbosityLevel | verbLevel = Teuchos::Describable::verbLevel_default |
||
) | const |
Print the object with some verbosity level to an FancyOStream object.