Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | List of all members
Ifpack2::LocalSparseTriangularSolver< MatrixType > Class Template Reference

"Preconditioner" that solves local sparse triangular systems. More...

#include <Ifpack2_LocalSparseTriangularSolver_decl.hpp>

Inheritance diagram for Ifpack2::LocalSparseTriangularSolver< MatrixType >:
Inheritance graph
[legend]

Public Types

typedef MatrixType::scalar_type scalar_type
 Type of the entries of the input matrix. More...
 
typedef
MatrixType::local_ordinal_type 
local_ordinal_type
 Type of the local indices of the input matrix. More...
 
typedef
MatrixType::global_ordinal_type 
global_ordinal_type
 Type of the global indices of the input matrix. More...
 
typedef MatrixType::node_type node_type
 Node type of the input matrix. More...
 
typedef MatrixType::mag_type magnitude_type
 Type of the absolute value (magnitude) of a scalar_type value. More...
 
typedef Tpetra::Map
< local_ordinal_type,
global_ordinal_type, node_type
map_type
 Specialization of Tpetra::Map used by this class. More...
 
typedef Tpetra::RowMatrix
< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type
row_matrix_type
 Specialization of Tpetra::RowMatrix used by this class. More...
 
typedef Tpetra::CrsMatrix
< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type
crs_matrix_type
 Specialization of Tpetra::CrsMatrix used by this class. 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

 LocalSparseTriangularSolver (const Teuchos::RCP< const row_matrix_type > &A)
 Constructor. More...
 
 LocalSparseTriangularSolver (const Teuchos::RCP< const row_matrix_type > &A, const Teuchos::RCP< Teuchos::FancyOStream > &out)
 Constructor that takes an optional debug output stream. More...
 
 LocalSparseTriangularSolver ()
 Constructor that takes no input matrix. More...
 
 LocalSparseTriangularSolver (const bool, const Teuchos::RCP< Teuchos::FancyOStream > &out)
 Constructor that takes no input matrix. More...
 
virtual ~LocalSparseTriangularSolver ()
 Destructor (virtual for memory safety). More...
 
void setParameters (const Teuchos::ParameterList &params)
 Set this object's parameters. More...
 
void initialize ()
 "Symbolic" phase of setup More...
 
bool isInitialized () const
 Return true if the preconditioner has been successfully initialized. More...
 
void compute ()
 "Numeric" phase of setup More...
 
bool isComputed () const
 Return true if compute() has been called. 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, and put the result in Y. More...
 
Teuchos::RCP< const map_typegetDomainMap () const
 The domain of this operator. More...
 
Teuchos::RCP< const map_typegetRangeMap () const
 The range of this operator. More...
 
void applyMat (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) const
 Apply the original input matrix. More...
 
Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm () const
 This operator's communicator. More...
 
Teuchos::RCP< const
row_matrix_type
getMatrix () const
 The original input matrix. More...
 
double getComputeFlops () const
 Return the number of flops in the computation phase. More...
 
double getApplyFlops () const
 Return the number of flops for the application of the preconditioner. More...
 
int getNumInitialize () const
 Return the number of calls to initialize(). More...
 
int getNumCompute () const
 Return the number of calls to compute(). More...
 
int getNumApply () const
 Return the number of calls to apply(). More...
 
double getInitializeTime () const
 Return the time spent in initialize(). More...
 
double getComputeTime () const
 Return the time spent in compute(). More...
 
double getApplyTime () const
 Return the time spent in apply(). More...
 
Implementation of Teuchos::Describable
std::string description () const
 A one-line description of this object. More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Print this object with given verbosity to the given output stream. More...
 
virtual void setMatrix (const Teuchos::RCP< const row_matrix_type > &A)
 Set this preconditioner's matrix. 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...
 
- Public Member Functions inherited from Ifpack2::Details::CanChangeMatrix< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > >
virtual void setMatrix (const Teuchos::RCP< const Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type > > &A)=0
 Set the new matrix. More...
 
virtual ~CanChangeMatrix ()
 Destructor. More...
 

Detailed Description

template<class MatrixType>
class Ifpack2::LocalSparseTriangularSolver< MatrixType >

"Preconditioner" that solves local sparse triangular systems.

Template Parameters
MatrixTypeSpecialization of Tpetra::RowMatrix.

This class solves local sparse triangular systems. "Local" means "per MPI process." The matrix itself may be distributed across multiple MPI processes, but this class works on each MPI process' part of the matrix, and the input and output multivectors, separately. (See this class' constructor for details.)

This effectively assumes that the global matrix is block diagonal. Triangular solves usually imply that these blocks are square. If a particular triangular solver knows how to deal with nonsquare blocks, though, this is allowed.

The implementation currently requires that the input Tpetra::RowMatrix actually be a Tpetra::CrsMatrix. This lets us optimize without necessarily copying data structures. We may relax this restriction in the future.

If you are writing a new Ifpack2 class that needs to solve local sparse triangular systems stored as Tpetra::CrsMatrix, use this class only.

Member Typedef Documentation

template<class MatrixType >
typedef MatrixType::scalar_type Ifpack2::LocalSparseTriangularSolver< MatrixType >::scalar_type

Type of the entries of the input matrix.

template<class MatrixType >
typedef MatrixType::local_ordinal_type Ifpack2::LocalSparseTriangularSolver< MatrixType >::local_ordinal_type

Type of the local indices of the input matrix.

template<class MatrixType >
typedef MatrixType::global_ordinal_type Ifpack2::LocalSparseTriangularSolver< MatrixType >::global_ordinal_type

Type of the global indices of the input matrix.

template<class MatrixType >
typedef MatrixType::node_type Ifpack2::LocalSparseTriangularSolver< MatrixType >::node_type

Node type of the input matrix.

template<class MatrixType >
typedef MatrixType::mag_type Ifpack2::LocalSparseTriangularSolver< MatrixType >::magnitude_type

Type of the absolute value (magnitude) of a scalar_type value.

template<class MatrixType >
typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> Ifpack2::LocalSparseTriangularSolver< MatrixType >::map_type

Specialization of Tpetra::Map used by this class.

template<class MatrixType >
typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::LocalSparseTriangularSolver< MatrixType >::row_matrix_type

Specialization of Tpetra::RowMatrix used by this class.

template<class MatrixType >
typedef Tpetra::CrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::LocalSparseTriangularSolver< MatrixType >::crs_matrix_type

Specialization of Tpetra::CrsMatrix used by this class.

Constructor & Destructor Documentation

template<class MatrixType >
Ifpack2::LocalSparseTriangularSolver< MatrixType >::LocalSparseTriangularSolver ( const Teuchos::RCP< const row_matrix_type > &  A)

Constructor.

Parameters
A[in] The input sparse matrix. Though its type is Tpetra::RowMatrix for consistency with other Ifpack2 solvers, this must be a Tpetra::CrsMatrix specialization.

The input matrix A may be distributed across multiple MPI processes. This class' apply() method will use A's Import object, if it exists, to Import the input MultiVector from the domain Map to the column Map. It will also use A's Export object, if it exists, to Export the output MultiVector from the row Map to the range Map. Thus, to avoid MPI communication and local permutations, construct A so that the row, column, range, and domain Maps are all identical.

On the other hand, you may encode local permutations in the matrix's Maps, and let Import and/or Export execute them for you.

The input matrix must have local properties corresponding to the way in which one wants to solve. ("Local" means "to each MPI process.") For example, if one wants to solve lower triangular systems with an implicit unit diagonal, the matrix A must have these properties. If the matrix does not know whether it has these properties and the user does not specify them, then this class is responsible for figuring out whether the matrix has those properties.

template<class MatrixType >
Ifpack2::LocalSparseTriangularSolver< MatrixType >::LocalSparseTriangularSolver ( const Teuchos::RCP< const row_matrix_type > &  A,
const Teuchos::RCP< Teuchos::FancyOStream > &  out 
)

Constructor that takes an optional debug output stream.

Parameters
A[in] The input sparse matrix. Though its type is Tpetra::RowMatrix for consistency with other Ifpack2 solvers, this must be a Tpetra::CrsMatrix specialization.
out[in/out] Optional debug output stream. If nonnull, this solver will print copious debug output to the stream.
template<class MatrixType >
Ifpack2::LocalSparseTriangularSolver< MatrixType >::LocalSparseTriangularSolver ( )

Constructor that takes no input matrix.

Call setMatrix to initialize the matrix.

template<class MatrixType >
Ifpack2::LocalSparseTriangularSolver< MatrixType >::LocalSparseTriangularSolver ( const bool  ,
const Teuchos::RCP< Teuchos::FancyOStream > &  out 
)

Constructor that takes no input matrix.

Call setMatrix to initialize the matrix.

Parameters
unused[in] Disambiguate the constructor so that the ctor taking A can cast A implicitly. This ctor is meant for debugging, so I prefer to make it a bit awkward than require careful explicit casting of A in the other ctor.
out[in/out] Optional debug output stream. If nonnull, this solver will print copious debug output to the stream.
template<class MatrixType >
Ifpack2::LocalSparseTriangularSolver< MatrixType >::~LocalSparseTriangularSolver ( )
virtual

Destructor (virtual for memory safety).

Member Function Documentation

template<class MatrixType >
void Ifpack2::LocalSparseTriangularSolver< MatrixType >::setParameters ( const Teuchos::ParameterList params)
virtual

Set this object's parameters.

Parameters
plist[in] List of parameters.
  • "trisolver: reverse U" (bool): reverse storage for upper triangular matrices to be more cache-efficient

If Trilinos_ENABLE_ShyLU_NodeHTS=TRUE, then these parameters are available:

  • "trisolver: type" (string): One of {"Internal" (default), "HTS"}.
  • "trisolver: block size" (int): The triangular matrix can usefully be thought of as being blocked int little blocks of this size. Default is 1.

Implements Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >.

template<class MatrixType >
void Ifpack2::LocalSparseTriangularSolver< MatrixType >::initialize ( )
virtual

"Symbolic" phase of setup

Call this before calling compute() or apply() if the matrix object itself changes, or if the matrix's graph structure may have changed.

Implements Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >.

template<class MatrixType >
bool Ifpack2::LocalSparseTriangularSolver< MatrixType >::isInitialized ( ) const
inlinevirtual
template<class MatrixType >
void Ifpack2::LocalSparseTriangularSolver< MatrixType >::compute ( )
virtual

"Numeric" phase of setup

Call this before calling apply() if the values in the matrix may have changed.

Implements Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >.

template<class MatrixType >
bool Ifpack2::LocalSparseTriangularSolver< MatrixType >::isComputed ( ) const
inlinevirtual
template<class MatrixType >
void Ifpack2::LocalSparseTriangularSolver< MatrixType >::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
virtual

Apply the preconditioner to X, and put the result in Y.

If this preconditioner is an operator M, this method computes Y := beta * Y + alpha * (M * X) .

Parameters
X[in] MultiVector to which to apply the preconditioner.
Y[in/out] On input: Initial guess, if applicable. On output: Result of applying the preconditioner. Y may alias X.
mode[in] Whether to apply the transpose (Teuchos::TRANS) or conjugate transpose (Teuchos::CONJ_TRANS). The default (Teuchos::NO_TRANS) is not to apply the transpose or conjugate transpose.
alpha[in] Scalar factor by which to multiply the result of applying this operator to X.
beta[in] Scalar factor for Y.

Implements Ifpack2::Preconditioner< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >.

template<class MatrixType >
Teuchos::RCP< const typename LocalSparseTriangularSolver< MatrixType >::map_type > Ifpack2::LocalSparseTriangularSolver< MatrixType >::getDomainMap ( ) const
virtual
template<class MatrixType >
Teuchos::RCP< const typename LocalSparseTriangularSolver< MatrixType >::map_type > Ifpack2::LocalSparseTriangularSolver< MatrixType >::getRangeMap ( ) const
virtual
template<class MatrixType >
void Ifpack2::LocalSparseTriangularSolver< MatrixType >::applyMat ( 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 
) const

Apply the original input matrix.

Parameters
X[in] MultiVector input.
Y[in/out] Result of applying Op(A) to X, where Op(A) is A (mode == Teuchos::NO_TRANS), the transpose of A (mode == Teuchos::TRANS), or the conjugate transpose of A (mode == Teuchos::CONJ_TRANS)
mode[in] See above.
template<class MatrixType >
Teuchos::RCP<const Teuchos::Comm<int> > Ifpack2::LocalSparseTriangularSolver< MatrixType >::getComm ( ) const

This operator's communicator.

template<class MatrixType >
Teuchos::RCP<const row_matrix_type> Ifpack2::LocalSparseTriangularSolver< MatrixType >::getMatrix ( ) const
inlinevirtual
template<class MatrixType >
double Ifpack2::LocalSparseTriangularSolver< MatrixType >::getComputeFlops ( ) const

Return the number of flops in the computation phase.

template<class MatrixType >
double Ifpack2::LocalSparseTriangularSolver< MatrixType >::getApplyFlops ( ) const

Return the number of flops for the application of the preconditioner.

template<class MatrixType >
int Ifpack2::LocalSparseTriangularSolver< MatrixType >::getNumInitialize ( ) const
virtual
template<class MatrixType >
int Ifpack2::LocalSparseTriangularSolver< MatrixType >::getNumCompute ( ) const
virtual
template<class MatrixType >
int Ifpack2::LocalSparseTriangularSolver< MatrixType >::getNumApply ( ) const
virtual
template<class MatrixType >
double Ifpack2::LocalSparseTriangularSolver< MatrixType >::getInitializeTime ( ) const
virtual
template<class MatrixType >
double Ifpack2::LocalSparseTriangularSolver< MatrixType >::getComputeTime ( ) const
virtual
template<class MatrixType >
double Ifpack2::LocalSparseTriangularSolver< MatrixType >::getApplyTime ( ) const
virtual
template<class MatrixType >
std::string Ifpack2::LocalSparseTriangularSolver< MatrixType >::description ( ) const

A one-line description of this object.

template<class MatrixType >
void Ifpack2::LocalSparseTriangularSolver< MatrixType >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const

Print this object with given verbosity to the given output stream.

Parameters
out[out] Output stream to which to print
verbLevel[in] Verbosity level

You may create a Teuchos::FancyOStream from any std::ostream. For example, to wrap std::cout in a FancyOStream, do this:

Teuchos::getFancyOStream (Teuchos::rcpFromRef (std::cout));

To wrap a new std::ostringstream in a FancyOStream, do this:

auto osPtr = Teuchos::rcp (new std::ostringstream ());
Teuchos::RCP<Teuchos::FancyOStream> out = Teuchos::getFancyOStream (osPtr);
// ... use out ...
// Print to std::cout whatever the std::ostringstream got.
std::cout << osPtr->str () << std::endl;
template<class MatrixType >
void Ifpack2::LocalSparseTriangularSolver< MatrixType >::setMatrix ( const Teuchos::RCP< const row_matrix_type > &  A)
virtual

Set this preconditioner's matrix.

After calling this method, you must call first initialize(), then compute(), before you may call apply().


The documentation for this class was generated from the following files: