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::IlukGraph< GraphType > Class Template Reference

Construct a level filled graph for use in computing an ILU(k) incomplete factorization. More...

#include <Ifpack2_IlukGraph.hpp>

Inheritance diagram for Ifpack2::IlukGraph< GraphType >:
Inheritance graph
[legend]

Public Types

typedef Tpetra::RowGraph
< local_ordinal_type,
global_ordinal_type, node_type > 
row_graph_type
 Tpetra::RowGraph specialization used by this class. More...
 
typedef Tpetra::CrsGraph
< local_ordinal_type,
global_ordinal_type, node_type > 
crs_graph_type
 Tpetra::CrsGraph specialization used by this class. More...
 

Public Member Functions

 IlukGraph (const Teuchos::RCP< const GraphType > &G, const int levelFill, const int levelOverlap, const double overalloc=2.)
 Constructor. More...
 
virtual ~IlukGraph ()
 IlukGraph Destructor. More...
 
void setParameters (const Teuchos::ParameterList &parameterlist)
 Set parameters. More...
 
void initialize ()
 Set up the graph structure of the L and U factors. More...
 
int getLevelFill () const
 The level of fill used to construct this graph. More...
 
int getLevelOverlap () const
 The level of overlap used to construct this graph. More...
 
Teuchos::RCP< crs_graph_typegetL_Graph () const
 Returns the graph of lower triangle of the ILU(k) graph as a Tpetra::CrsGraph. More...
 
Teuchos::RCP< crs_graph_typegetU_Graph () const
 Returns the graph of upper triangle of the ILU(k) graph as a Tpetra::CrsGraph. More...
 
Teuchos::RCP< const
crs_graph_type
getOverlapGraph () const
 Returns the the overlapped graph. More...
 
size_t getNumGlobalDiagonals () const
 Returns the global number of diagonals in the ILU(k) graph. More...
 

Detailed Description

template<class GraphType>
class Ifpack2::IlukGraph< GraphType >

Construct a level filled graph for use in computing an ILU(k) incomplete factorization.

Template Parameters
GraphTypeA specialization of Tpetra::CrsGraph (tested) or Tpetra::RowGraph (not tested).
Warning
This class is an implementation detail of RILUK. Its interface may change or it may go away at any time.

Ifpack2::IlukGraph enables the construction of matrix graphs using level-fill algorithms. The only function required for construction is a getGlobalRowView capability, i.e., the graph that is passed in to the constructor must implement the Tpetra::RowGraph interface defined in Tpetra_RowGraph.hpp.

Constructing IlukGraph objects

Constructing an IlukGraph is a two step process. First, call the constructor, passing in a CrsGraph object and an integer indicating the desired level of fill. Then, call the initialize() method to complete the process. This naturally matches the three-stage initialization of Preconditioner objects, and in particular of RILUK.

It is worth noting that an IlukGraph object creates two Tpetra::CrsGraph objects containing L and U, the graphs for the lower and upper triangular parts of the ILU(k) graph. Thus, it is possible to manually insert and delete graph entries in L and U via the Tpetra::CrsGraph InsertIndices and RemoveIndices functions. However, in this case FillComplete must be called before the graph is used for subsequent operations.

Member Typedef Documentation

template<class GraphType >
typedef Tpetra::RowGraph<local_ordinal_type, global_ordinal_type, node_type> Ifpack2::IlukGraph< GraphType >::row_graph_type

Tpetra::RowGraph specialization used by this class.

template<class GraphType >
typedef Tpetra::CrsGraph<local_ordinal_type, global_ordinal_type, node_type> Ifpack2::IlukGraph< GraphType >::crs_graph_type

Tpetra::CrsGraph specialization used by this class.

Constructor & Destructor Documentation

template<class GraphType >
Ifpack2::IlukGraph< GraphType >::IlukGraph ( const Teuchos::RCP< const GraphType > &  G,
const int  levelFill,
const int  levelOverlap,
const double  overalloc = 2. 
)

Constructor.

Create a IlukGraph object using the input graph and specified level of fill.

Parameters
G[in] An existing graph.
levelFill[in] The level of fill to compute; the k of ILU(k).
levelOverlap[in] The level of overlap between subdomains.
overalloc[in] The estimated number of nonzeros per row in the resulting matrices is (maxNodeNumRowEntries of the input * overalloc^levelFill). Must be greater than 1. Smaller values are more conservative with memory, but may require recomputation if the estimate is too low. Default value is two.
Note
Actual construction occurs in initialize().
template<class GraphType >
Ifpack2::IlukGraph< GraphType >::~IlukGraph ( )
virtual

IlukGraph Destructor.

Member Function Documentation

template<class GraphType >
void Ifpack2::IlukGraph< GraphType >::setParameters ( const Teuchos::ParameterList parameterlist)

Set parameters.

This method recognizes two parameter names: Level_fill and Level_overlap. Both are case insensitive, and in both cases the parameter must have type int.

template<class GraphType >
void Ifpack2::IlukGraph< GraphType >::initialize ( )

Set up the graph structure of the L and U factors.

This method is called "initialize" by analogy with Preconditioner, where initialize() computes the symbolic (incomplete) factorization, and compute() computes the corresponding numeric factorization. IlukGraph is just a graph, so it can only compute a symbolic factorization (i.e., the graph structure of the factorization). Hence, it implements initialize(), but not compute(). RILUK calls IlukGraph's initialize() method in its own initialize() method, as one would expect.

template<class GraphType >
int Ifpack2::IlukGraph< GraphType >::getLevelFill ( ) const
inline

The level of fill used to construct this graph.

template<class GraphType >
int Ifpack2::IlukGraph< GraphType >::getLevelOverlap ( ) const
inline

The level of overlap used to construct this graph.

template<class GraphType >
Teuchos::RCP<crs_graph_type> Ifpack2::IlukGraph< GraphType >::getL_Graph ( ) const
inline

Returns the graph of lower triangle of the ILU(k) graph as a Tpetra::CrsGraph.

template<class GraphType >
Teuchos::RCP<crs_graph_type> Ifpack2::IlukGraph< GraphType >::getU_Graph ( ) const
inline

Returns the graph of upper triangle of the ILU(k) graph as a Tpetra::CrsGraph.

template<class GraphType >
Teuchos::RCP<const crs_graph_type> Ifpack2::IlukGraph< GraphType >::getOverlapGraph ( ) const
inline

Returns the the overlapped graph.

template<class GraphType >
size_t Ifpack2::IlukGraph< GraphType >::getNumGlobalDiagonals ( ) const
inline

Returns the global number of diagonals in the ILU(k) graph.


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