Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Types | Public Member Functions | Protected Member Functions | List of all members
Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Referenceabstract

The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic) More...

#include <Ifpack2_Details_FastILU_Base_decl.hpp>

Inheritance diagram for Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >:
Inheritance graph
[legend]

Public Types

typedef Node::device_type device_type
 Kokkos device type. More...
 
typedef
device_type::execution_space 
execution_space
 Kokkos execution space. More...
 
typedef Tpetra::MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node >
::impl_scalar_type 
ImplScalar
 Kokkos scalar type. More...
 
typedef Tpetra::RowMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > 
TRowMatrix
 Tpetra row matrix. More...
 
typedef Tpetra::CrsMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > 
TCrsMatrix
 Tpetra CRS matrix. More...
 
typedef Tpetra::MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > 
TMultiVec
 Tpetra multivector. More...
 
typedef
KokkosSparse::CrsMatrix
< Scalar, LocalOrdinal,
execution_space
KCrsMatrix
 Kokkos CRS matrix. More...
 
typedef Kokkos::View
< LocalOrdinal
*, execution_space
OrdinalArray
 Array of LocalOrdinal on device. More...
 
typedef Kokkos::View
< LocalOrdinal
*, execution_space >
::HostMirror 
OrdinalArrayHost
 Array of LocalOrdinal on host. More...
 
typedef Kokkos::View
< ImplScalar
*, execution_space
ImplScalarArray
 Array of Scalar on device. More...
 
- Public Types inherited from Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >
typedef Teuchos::ScalarTraits
< Scalar >::magnitudeType 
magnitude_type
 The type of the magnitude (absolute value) of a matrix entry. More...
 

Public Member Functions

 FastILU_Base (Teuchos::RCP< const TRowMatrix > mat_)
 Constructor. More...
 
Teuchos::RCP< const
Tpetra::Map< LocalOrdinal,
GlobalOrdinal, Node > > 
getDomainMap () const
 Get the domain map of the matrix. More...
 
Teuchos::RCP< const
Tpetra::Map< LocalOrdinal,
GlobalOrdinal, Node > > 
getRangeMap () const
 Get the range map of the matrix. More...
 
void apply (const TMultiVec &X, TMultiVec &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const
 Apply the preconditioner. More...
 
void setParameters (const Teuchos::ParameterList &List)
 Validate parameters, and set defaults when parameters are not provided. More...
 
void initialize ()
 Initialize the preconditioner. More...
 
bool isInitialized () const
 Whether initialize() has been called since the last time the matrix's structure was changed. More...
 
void compute ()
 Compute the preconditioner. More...
 
bool isComputed () const
 Whether compute() has been called since the last time the matrix's values or structure were changed. More...
 
Teuchos::RCP< const TRowMatrixgetMatrix () const
 Get the current matrix. More...
 
int getNumInitialize () const
 Get the number of times initialize() was called. More...
 
int getNumCompute () const
 Get the number of times compute() was called. More...
 
int getNumApply () const
 Get the number of times apply() was called. More...
 
double getInitializeTime () const
 Get the time spent in the last initialize() call. More...
 
double getComputeTime () const
 Get the time spent in the last compute() call. More...
 
double getApplyTime () const
 Get the time spent in the last apply() call. More...
 
double getCopyTime () const
 Get the time spent deep copying local 3-array CRS out of the matrix. More...
 
virtual int getSweeps () const =0
 Get the "sweeps" parameter. More...
 
virtual std::string getSpTrsvType () const =0
 Get the name of triangular solve algorithm. More...
 
virtual int getNTrisol () const =0
 Get the "triangular solve iterations" parameter. More...
 
virtual void checkLocalILU () const
 Verify and print debug information about the underlying ILU preconditioner (only supported if this is an Ifpack2::Details::Filu) More...
 
virtual void checkLocalIC () const
 Verify and print debug information about the underlying IC preconditioner. More...
 
std::string description () const
 Return a brief description of the preconditioner, in YAML format. More...
 
void setMatrix (const Teuchos::RCP< const TRowMatrix > &A)
 
- Public Member Functions inherited from Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >
virtual ~Preconditioner ()
 Destructor. More...
 
virtual void setZeroStartingSolution (bool zeroStartingSolution)
 Set this preconditioner's parameters. More...
 
- Public Member Functions inherited from Ifpack2::Details::CanChangeMatrix< Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > >
virtual void setMatrix (const Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A)=0
 Set the new matrix. More...
 
virtual ~CanChangeMatrix ()
 Destructor. More...
 

Protected Member Functions

virtual void initLocalPrec ()=0
 Construct the underlying preconditioner (localPrec_) using given params and then call localPrec_->initialize() More...
 
virtual void computeLocalPrec ()=0
 Get values array from the matrix and then call compute() on the underlying preconditioner. More...
 
virtual void applyLocalPrec (ImplScalarArray x, ImplScalarArray y) const =0
 Apply the local preconditioner with 1-D views of the local parts of X and Y (one vector only) More...
 
virtual std::string getName () const =0
 Get the name of the underlying preconditioner ("Filu", "Fildl" or "Fic") More...
 

Detailed Description

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
class Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >

The base class of the Ifpack2 FastILU wrappers (Filu, Fildl and Fic)

Member Typedef Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Node::device_type Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::device_type

Kokkos device type.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef device_type::execution_space Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::execution_space

Kokkos execution space.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Tpetra::MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node>::impl_scalar_type Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ImplScalar

Kokkos scalar type.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Tpetra::RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::TRowMatrix

Tpetra row matrix.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::TCrsMatrix

Tpetra CRS matrix.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::TMultiVec

Tpetra multivector.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef KokkosSparse::CrsMatrix<Scalar, LocalOrdinal, execution_space> Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::KCrsMatrix

Kokkos CRS matrix.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Kokkos::View<LocalOrdinal *, execution_space> Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::OrdinalArray

Array of LocalOrdinal on device.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Kokkos::View<LocalOrdinal *, execution_space>::HostMirror Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::OrdinalArrayHost

Array of LocalOrdinal on host.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
typedef Kokkos::View< ImplScalar *, execution_space> Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ImplScalarArray

Array of Scalar on device.

Constructor & Destructor Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::FastILU_Base ( Teuchos::RCP< const TRowMatrix mat_)

Constructor.

Member Function Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getDomainMap ( ) const
virtual

Get the domain map of the matrix.

Implements Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getRangeMap ( ) const
virtual

Get the range map of the matrix.

Implements Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::apply ( const TMultiVec X,
TMultiVec Y,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
Scalar  alpha = Teuchos::ScalarTraits<Scalar>::one(),
Scalar  beta = Teuchos::ScalarTraits<Scalar>::zero() 
) const
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setParameters ( const Teuchos::ParameterList List)
virtual

Validate parameters, and set defaults when parameters are not provided.

Available parameters:

Parameter Name Parameter Type Default Description
"sweeps" int 5 The number of iterations of the preconditioner during apply()
"triangular solve iterations" int 1 The number of iterations of the block Jacobi triangular solver during apply()
"level" int 0 The level of fill
"damping factor" double 0.5 The damping factor used during apply() – must be between 0 (inclusive) and 1 (exclusive)
"shift" double 0 The Manteuffel shift parameter
"guess" bool true Whether to create multiple preconditioners at successively lower levels of fill to create the initial guess
"block size" int 1 The block size for the block Jacobi iterations

Implements Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::initialize ( )
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isInitialized ( ) const
virtual

Whether initialize() has been called since the last time the matrix's structure was changed.

Implements Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::compute ( )
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
bool Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::isComputed ( ) const
virtual

Whether compute() has been called since the last time the matrix's values or structure were changed.

Implements Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Tpetra::RowMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getMatrix ( ) const
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumInitialize ( ) const
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumCompute ( ) const
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
int Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNumApply ( ) const
virtual

Get the number of times apply() was called.

Implements Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
double Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getInitializeTime ( ) const
virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
double Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getComputeTime ( ) const
virtual

Get the time spent in the last compute() call.

Implements Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
double Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getApplyTime ( ) const
virtual

Get the time spent in the last apply() call.

Implements Ifpack2::Preconditioner< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
double Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getCopyTime ( ) const

Get the time spent deep copying local 3-array CRS out of the matrix.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
virtual int Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getSweeps ( ) const
pure virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
virtual std::string Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getSpTrsvType ( ) const
pure virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
virtual int Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getNTrisol ( ) const
pure virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::checkLocalILU ( ) const
virtual

Verify and print debug information about the underlying ILU preconditioner (only supported if this is an Ifpack2::Details::Filu)

Reimplemented in Ifpack2::Details::Filu< Scalar, LocalOrdinal, GlobalOrdinal, Node, BlockCrsEnabled >.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::checkLocalIC ( ) const
virtual
template<typename Scalar , typename LocalOrdinal , typename GlobalOrdinal , typename Node >
std::string Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::description ( ) const

Return a brief description of the preconditioner, in YAML format.

template<typename Scalar , typename LocalOrdinal , typename GlobalOrdinal , typename Node >
void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::setMatrix ( const Teuchos::RCP< const TRowMatrix > &  A)

Provide a new matrix If the A's graph is different from that of the existing matrix, initialize() and compute() will need to be called again before the next apply(). If A's values are different than those of the existing matrix, compute() will need to be called again before the next apply().

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
virtual void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::initLocalPrec ( )
protectedpure virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
virtual void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::computeLocalPrec ( )
protectedpure virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
virtual void Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::applyLocalPrec ( ImplScalarArray  x,
ImplScalarArray  y 
) const
protectedpure virtual
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
virtual std::string Ifpack2::Details::FastILU_Base< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getName ( ) const
protectedpure virtual

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