Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Static Public Member Functions | List of all members
Ifpack2::Factory Class Reference

"Factory" for creating Ifpack2 preconditioners. More...

#include <Ifpack2_Factory_decl.hpp>

Static Public Member Functions

template<class MatrixType >
static Teuchos::RCP
< Preconditioner< typename
MatrixType::scalar_type,
typename
MatrixType::local_ordinal_type,
typename
MatrixType::global_ordinal_type,
typename MatrixType::node_type > > 
create (const std::string &precType, const Teuchos::RCP< const MatrixType > &matrix)
 Create an instance of Ifpack2_Preconditioner given the string name of the preconditioner type. More...
 
template<class MatrixType >
static Teuchos::RCP
< Preconditioner< typename
MatrixType::scalar_type,
typename
MatrixType::local_ordinal_type,
typename
MatrixType::global_ordinal_type,
typename MatrixType::node_type > > 
create (const std::string &precType, const Teuchos::RCP< const MatrixType > &matrix, const int overlap)
 Create an instance of Ifpack2_Preconditioner given the string name of the preconditioner type. More...
 
template<class InputMatrixType , class OutputMatrixType >
static Teuchos::RCP
< Preconditioner< typename
OutputMatrixType::scalar_type,
typename
OutputMatrixType::local_ordinal_type,
typename
OutputMatrixType::global_ordinal_type,
typename
OutputMatrixType::node_type > > 
clone (const Teuchos::RCP< Preconditioner< typename InputMatrixType::scalar_type, typename InputMatrixType::local_ordinal_type, typename InputMatrixType::global_ordinal_type, typename InputMatrixType::node_type > > &prec, const Teuchos::RCP< const OutputMatrixType > &matrix, const Teuchos::ParameterList &params=Teuchos::ParameterList())
 Clones a preconditioner for a different node type from an Ifpack2 RILUK or Chebyshev preconditioner. More...
 

Detailed Description

"Factory" for creating Ifpack2 preconditioners.

This class' create() method lets users create an instance of any Ifpack2 preconditioner.

The create() method has three arguments:

The first argument can assume the following values:

The following fragment of code shows the basic usage of this class.

#include "Ifpack2_Factory.hpp"
// ...
typedef double Scalar;
typedef Tpetra::CrsMatrix<Scalar> crs_matrix_type;
// ...
RCP<crs_matrix_type> A;
// ... fill A and call fillComplete() ...
// Use ILUT (incomplete LU with thresholding) on each process.
const std::string precType = "ILUT";
RCP<prec_type> prec = factory.create (precType, A);
ParameterList params;
params.set ("fact: ilut level-of-fill", 5.0); // ILUT(fill=5, drop=0)
prec->setParameters (params);
prec->initialize ();
prec->compute ();
// now prec can be used as a preconditioner

Member Function Documentation

template<class MatrixType >
static Teuchos::RCP<Preconditioner<typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type> > Ifpack2::Factory::create ( const std::string &  precType,
const Teuchos::RCP< const MatrixType > &  matrix 
)
inlinestatic

Create an instance of Ifpack2_Preconditioner given the string name of the preconditioner type.

Parameters
precType[in] Name of preconditioner type to be created.
matrix[in] Matrix used to define the preconditioner

Throw an exception if the preconditioner with that input name does not exist. Otherwise, return a newly created preconditioner object.

template<class MatrixType >
static Teuchos::RCP<Preconditioner<typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type> > Ifpack2::Factory::create ( const std::string &  precType,
const Teuchos::RCP< const MatrixType > &  matrix,
const int  overlap 
)
inlinestatic

Create an instance of Ifpack2_Preconditioner given the string name of the preconditioner type.

Warning
This version of the constructor is DEPRECATED, because the single-argument version suffices; users may specify the overlap level via the "schwarz: overlap level" parameter.
Parameters
precType[in] Name of preconditioner type to be created.
matrix[in] Matrix used to define the preconditioner
overlap(in) Specified overlap; defaults to 0.

Throw an exception if the preconditioner with that input name does not exist. Otherwise, return a newly created preconditioner object.

template<class InputMatrixType , class OutputMatrixType >
static Teuchos::RCP<Preconditioner<typename OutputMatrixType::scalar_type, typename OutputMatrixType::local_ordinal_type, typename OutputMatrixType::global_ordinal_type, typename OutputMatrixType::node_type> > Ifpack2::Factory::clone ( const Teuchos::RCP< Preconditioner< typename InputMatrixType::scalar_type, typename InputMatrixType::local_ordinal_type, typename InputMatrixType::global_ordinal_type, typename InputMatrixType::node_type > > &  prec,
const Teuchos::RCP< const OutputMatrixType > &  matrix,
const Teuchos::ParameterList params = Teuchos::ParameterList () 
)
inlinestatic

Clones a preconditioner for a different node type from an Ifpack2 RILUK or Chebyshev preconditioner.


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