Ifpack2 Templated Preconditioning Package Version 1.0
Public Member Functions
Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node > Class Template Reference

A class for constructing level filled graphs for use with ILU(k) class preconditioners. More...

#include <Ifpack2_IlukGraph.hpp>

Inheritance diagram for Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >:
Inheritance graph
[legend]

List of all members.

Public Member Functions

 IlukGraph (const Teuchos::RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &Graph_in, int LevelFill_in, int LevelOverlap_in)
 IlukGraph constuctor.
 IlukGraph (const IlukGraph< LocalOrdinal, GlobalOrdinal, Node > &Graph_in)
 Copy constructor.
virtual ~IlukGraph ()
 IlukGraph Destructor.
void setParameters (const Teuchos::ParameterList &parameterlist)
 Set parameters using Teuchos::ParameterList object.
void constructFilledGraph ()
 This method performs the actual construction of the graph.
void constructOverlapGraph ()
 Does the actual construction of the overlap matrix graph.
int getLevelFill () const
 Returns the level of fill used to construct this graph.
int getLevelOverlap () const
 Returns the level of overlap used to construct this graph.
const Teuchos::RCP
< Tpetra::CrsGraph
< LocalOrdinal, GlobalOrdinal,
Node > > & 
getL_Graph () const
 Returns the graph of lower triangle of the ILU(k) graph as a Tpetra::CrsGraph.
const Teuchos::RCP
< Tpetra::CrsGraph
< LocalOrdinal, GlobalOrdinal,
Node > > & 
getU_Graph () const
 Returns the graph of upper triangle of the ILU(k) graph as a Tpetra::CrsGraph.
const Teuchos::RCP
< Tpetra::CrsGraph
< LocalOrdinal, GlobalOrdinal,
Node > > & 
getOverlapGraph () const
 Returns the the overlapped graph.
size_t getNumGlobalDiagonals () const
 Returns the global number of diagonals in the ILU(k) graph.

Detailed Description

template<class LocalOrdinal, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
class Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >

A class for constructing level filled graphs for use with ILU(k) class preconditioners.

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 RowGraph interface defined in Tpetra_RowGraph.hpp

Constructing IlukGraph objects

Constructing IlukGraph objects is usually a two step process of passing in a CrsGraph object and an integer indicating the desired level of fill and then calling the ConstructFilledGraph function to complete the process. This allows warning error codes to be returned to the calling routine.

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.


Constructor & Destructor Documentation

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::IlukGraph ( const Teuchos::RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &  Graph_in,
int  LevelFill_in,
int  LevelOverlap_in 
)

IlukGraph constuctor.

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

Parameters:
InGraph_in - An existing CrsGraph. This object must implement the CrsGraph functions that provide graph dimension and pattern information.
InLevelFill_in - The level of fill to compute via ILU(k) algorithm.
InLevelOverlap_in - The level of overlap between subdomains.
Warning:
Actual construction occurs in constructFilledGraph.
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::IlukGraph ( const IlukGraph< LocalOrdinal, GlobalOrdinal, Node > &  Graph_in)

Copy constructor.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::~IlukGraph ( ) [virtual]

IlukGraph Destructor.


Member Function Documentation

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::setParameters ( const Teuchos::ParameterList parameterlist)

Set parameters using Teuchos::ParameterList object.

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

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::constructFilledGraph ( )

This method performs the actual construction of the graph.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::constructOverlapGraph ( )

Does the actual construction of the overlap matrix graph.

template<class LocalOrdinal, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
int Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::getLevelFill ( ) const [inline]

Returns the level of fill used to construct this graph.

template<class LocalOrdinal, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
int Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::getLevelOverlap ( ) const [inline]

Returns the level of overlap used to construct this graph.

template<class LocalOrdinal, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node> >& Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::getL_Graph ( ) const [inline]

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

template<class LocalOrdinal, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node> >& Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::getU_Graph ( ) const [inline]

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

template<class LocalOrdinal, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
const Teuchos::RCP<Tpetra::CrsGraph<LocalOrdinal,GlobalOrdinal,Node> >& Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::getOverlapGraph ( ) const [inline]

Returns the the overlapped graph.

template<class LocalOrdinal, class GlobalOrdinal = LocalOrdinal, class Node = Kokkos::DefaultNode::DefaultNodeType>
size_t Ifpack2::IlukGraph< LocalOrdinal, GlobalOrdinal, Node >::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:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends