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

Overlapping Schwarz where redundant patches are not stored explicitly. More...

#include <Ifpack2_DatabaseSchwarz_decl.hpp>

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

Public Types

Typedefs
typedef MatrixType matrix_type
 The template parameter of this class. More...
 
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::device_type 
device_type
 The Kokkos::Device specialization used by 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
 The Tpetra::RowMatrix specialization matching MatrixType. More...
 
typedef Tpetra::Map
< local_ordinal_type,
global_ordinal_type, node_type
map_type
 The Tpetra::Map specialization matching MatrixType. More...
 
typedef Tpetra::Vector
< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type
vector_type
 The Tpetra::Vector specialization matching MatrixType. 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

 DatabaseSchwarz (const Teuchos::RCP< const row_matrix_type > &A)
 Constructor. More...
 
 DatabaseSchwarz (const Teuchos::RCP< const row_matrix_type > &A, Teuchos::ParameterList &params)
 Constructor. More...
 
virtual ~DatabaseSchwarz ()
 Destructor. More...
 
Preconditioner computation methods
void setParameters (const Teuchos::ParameterList &params)
 Set (or reset) parameters. More...
 
bool supportsZeroStartingSolution ()
 
void setZeroStartingSolution (bool zeroStartingSolution)
 Set this preconditioner's parameters. More...
 
void initialize ()
 Initialize the preconditioner. More...
 
bool isInitialized () const
 
void compute ()
 (Re)compute the left scaling, and (if applicable) estimate max and min eigenvalues of D_inv * A. More...
 
bool isComputed () const
 Whether compute() has been called at least once. More...
 
Implementation of Ifpack2::Details::CanChangeMatrix
virtual void setMatrix (const Teuchos::RCP< const row_matrix_type > &A)
 Change the matrix to be preconditioned. 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, returning the result in Y. Y = alpha*Op(A)*X + beta*Y. More...
 
Teuchos::RCP< const map_typegetDomainMap () const
 The Tpetra::Map representing the domain of this operator. More...
 
Teuchos::RCP< const map_typegetRangeMap () const
 The Tpetra::Map representing the range of this operator. More...
 
bool hasTransposeApply () const
 Whether it's possible to apply the transpose 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
 Compute Y = Op(A)*X, where Op(A) is either A, \(A^T\), or \(A^H\). More...
 
Implementation of Teuchos::Describable
std::string description () const
 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 a Teuchos::FancyOStream. 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...
 

Attribute accessor methods

Teuchos::RCP< const
row_matrix_type
A_
 The matrix for which this is a preconditioner. More...
 
Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm () const
 The communicator over which the matrix is distributed. More...
 
Teuchos::RCP< const
row_matrix_type
getMatrix () const
 The matrix for which this is a preconditioner. More...
 
Teuchos::RCP< const
Tpetra::CrsMatrix< scalar_type,
local_ordinal_type,
global_ordinal_type, node_type > > 
getCrsMatrix () const
 Attempt to return the matrix A as a Tpetra::CrsMatrix. More...
 
double getComputeFlops () const
 The total number of floating-point operations taken by all calls to compute(). More...
 
double getApplyFlops () const
 The total number of floating-point operations taken by all calls to apply(). More...
 
int getNumInitialize () const
 The total number of successful calls to initialize(). More...
 
int getNumCompute () const
 The total number of successful calls to compute(). More...
 
int getNumApply () const
 The total number of successful calls to apply(). More...
 
double getInitializeTime () const
 The total time spent in all calls to initialize(). More...
 
double getComputeTime () const
 The total time spent in all calls to compute(). More...
 
double getApplyTime () const
 The total time spent in all calls to apply(). More...
 
size_t getNodeSmootherComplexity () const
 Get a rough estimate of cost per iteration. More...
 

Detailed Description

template<class MatrixType>
class Ifpack2::DatabaseSchwarz< MatrixType >

Overlapping Schwarz where redundant patches are not stored explicitly.

Template Parameters
MatrixTypeThe type of matrix to use.

This class implements an overlapping Schwarz assuming it has already been provided overlapping data containers. Additionally, it is assumed there is at least one row with PatchSize_ nonzeros for each patch. This means that a matrix obtained from piecewise linear FEMs in 2d, for example, will not have rows with 4 nonzeros for each patch. However, piecewise quadratic FEMs in 2d will have one row with 9 nonzeros for each cell in the mesh, meaning all the patches will be detected. For a matrix corresponding to a Taylor-Hood discretization in 2d (quadratic velocities and linear pressures), the patch size would be 22.

The general algorithm proceeds as follows:

  1. The rows of A are analyzed sequentially for any row with num_entries == PatchSize_
  2. If the current row corresponds to a DOF that has already been "visited", it is skipped
  3. Then, all nonzero indices belonging to the row are marked as "visited"
  4. The local patch matrix is formed and compared to a database of previous patch matrices. If any patch matrix in the database has an l1 distance to the current patch matrix less than tol, the current patch is not matrix is not stored, but an index pointing it to a replacement is stored instead
  5. Finally, if this patch matrix is not sufficiently close to one that has already been seen, it is added to the database
  6. The compute phase then inverts only the patch matcies in the database, and the apply phase loops over all patches and applies the inverse of each appropriate patch matrix

In general, there is a noticeable speedup when using this method compared to a typical method. This speedup may be further improved by using more advanced linear algebra interfaces such as batched Kokkos solves instead of the current LAPACK approach.

Member Typedef Documentation

template<class MatrixType>
typedef MatrixType Ifpack2::DatabaseSchwarz< MatrixType >::matrix_type

The template parameter of this class.

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

The type of the entries of the input MatrixType.

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

The type of local indices in the input MatrixType.

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

The type of global indices in the input MatrixType.

template<class MatrixType>
typedef MatrixType::node_type::device_type Ifpack2::DatabaseSchwarz< MatrixType >::device_type

The Kokkos::Device specialization used by the input MatrixType.

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

The Node type used by the input MatrixType.

template<class MatrixType>
typedef Teuchos::ScalarTraits<scalar_type>::magnitudeType Ifpack2::DatabaseSchwarz< MatrixType >::magnitude_type

The type of the magnitude (absolute value) of a matrix entry.

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

The Tpetra::RowMatrix specialization matching MatrixType.

MatrixType must be a Tpetra::RowMatrix specialization. This typedef will always be a Tpetra::RowMatrix specialization.

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

The Tpetra::Map specialization matching MatrixType.

template<class MatrixType>
typedef Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> Ifpack2::DatabaseSchwarz< MatrixType >::vector_type

The Tpetra::Vector specialization matching MatrixType.

If you wish to supply setParameters() a precomputed vector of diagonal entries of the matrix, use a pointer to an object of this type.

Constructor & Destructor Documentation

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

Constructor.

Parameters
[in]AThe sparse matrix to which to apply DatabaseSchwarz iteration. The matrix A must be square, and its domain Map and range Map must be the same. The latter means that the vectors x and y in the sparse matrix-vector product y = A*x must both have the same distribution over process(es).
template<class MatrixType >
Ifpack2::DatabaseSchwarz< MatrixType >::DatabaseSchwarz ( const Teuchos::RCP< const row_matrix_type > &  A,
Teuchos::ParameterList params 
)

Constructor.

Parameters
[in]AThe sparse matrix to which to apply DatabaseSchwarz iteration. The matrix A must be square, and its domain Map and range Map must be the same. The latter means that the vectors x and y in the sparse matrix-vector product y = A*x must both have the same distribution over process(es).
[in]paramsThe parameterlist containing settings for the object, such as the patch size to search for.
template<class MatrixType >
Ifpack2::DatabaseSchwarz< MatrixType >::~DatabaseSchwarz ( )
virtual

Destructor.

Member Function Documentation

template<class MatrixType >
void Ifpack2::DatabaseSchwarz< MatrixType >::setParameters ( const Teuchos::ParameterList params)
virtual
template<class MatrixType >
void Ifpack2::DatabaseSchwarz< MatrixType >::setZeroStartingSolution ( bool  zeroStartingSolution)
virtual
template<class MatrixType >
void Ifpack2::DatabaseSchwarz< MatrixType >::initialize ( )
virtual

Initialize the preconditioner.

The compute() method will call initialize() automatically if it has not yet been called, so you do not normally need to call this. However, it is correct to call initialize() yourself, and compute() will not call it again if it already has been called.

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

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

(Re)compute the left scaling, and (if applicable) estimate max and min eigenvalues of D_inv * A.

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

template<class MatrixType>
bool Ifpack2::DatabaseSchwarz< MatrixType >::isComputed ( ) const
inlinevirtual
template<class MatrixType >
void Ifpack2::DatabaseSchwarz< MatrixType >::setMatrix ( const Teuchos::RCP< const row_matrix_type > &  A)
virtual

Change the matrix to be preconditioned.

template<class MatrixType >
void Ifpack2::DatabaseSchwarz< 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, returning the result in Y. Y = alpha*Op(A)*X + beta*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 DatabaseSchwarz< MatrixType >::map_type > Ifpack2::DatabaseSchwarz< MatrixType >::getDomainMap ( ) const
virtual
template<class MatrixType >
Teuchos::RCP< const typename DatabaseSchwarz< MatrixType >::map_type > Ifpack2::DatabaseSchwarz< MatrixType >::getRangeMap ( ) const
virtual
template<class MatrixType >
bool Ifpack2::DatabaseSchwarz< MatrixType >::hasTransposeApply ( ) const

Whether it's possible to apply the transpose of this operator.

template<class MatrixType >
void Ifpack2::DatabaseSchwarz< 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

Compute Y = Op(A)*X, where Op(A) is either A, \(A^T\), or \(A^H\).

template<class MatrixType >
Teuchos::RCP< const Teuchos::Comm< int > > Ifpack2::DatabaseSchwarz< MatrixType >::getComm ( ) const

The communicator over which the matrix is distributed.

template<class MatrixType >
Teuchos::RCP< const typename DatabaseSchwarz< MatrixType >::row_matrix_type > Ifpack2::DatabaseSchwarz< MatrixType >::getMatrix ( ) const
virtual
template<class MatrixType >
Teuchos::RCP< const Tpetra::CrsMatrix< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > Ifpack2::DatabaseSchwarz< MatrixType >::getCrsMatrix ( ) const

Attempt to return the matrix A as a Tpetra::CrsMatrix.

This class does not require that A be a Tpetra::CrsMatrix. If it is NOT, this method will return Teuchos::null.

template<class MatrixType >
double Ifpack2::DatabaseSchwarz< MatrixType >::getComputeFlops ( ) const

The total number of floating-point operations taken by all calls to compute().

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

The total number of floating-point operations taken by all calls to apply().

template<class MatrixType >
int Ifpack2::DatabaseSchwarz< MatrixType >::getNumInitialize ( ) const
virtual
template<class MatrixType >
int Ifpack2::DatabaseSchwarz< MatrixType >::getNumCompute ( ) const
virtual
template<class MatrixType >
int Ifpack2::DatabaseSchwarz< MatrixType >::getNumApply ( ) const
virtual
template<class MatrixType >
double Ifpack2::DatabaseSchwarz< MatrixType >::getInitializeTime ( ) const
virtual
template<class MatrixType >
double Ifpack2::DatabaseSchwarz< MatrixType >::getComputeTime ( ) const
virtual
template<class MatrixType >
double Ifpack2::DatabaseSchwarz< MatrixType >::getApplyTime ( ) const
virtual
template<class MatrixType >
size_t Ifpack2::DatabaseSchwarz< MatrixType >::getNodeSmootherComplexity ( ) const

Get a rough estimate of cost per iteration.

template<class MatrixType >
std::string Ifpack2::DatabaseSchwarz< MatrixType >::description ( ) const

A simple one-line description of this object.

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

Print the object with some verbosity level to a Teuchos::FancyOStream.

Member Data Documentation

template<class MatrixType>
Teuchos::RCP<const row_matrix_type> Ifpack2::DatabaseSchwarz< MatrixType >::A_

The matrix for which this is a preconditioner.


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