Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class preconditioners. More...
#include <Ifpack_IlukGraph.h>
Public Member Functions | |
Ifpack_IlukGraph (const Epetra_CrsGraph &Graph_in, int LevelFill_in, int LevelOverlap_in) | |
Ifpack_IlukGraph constuctor. More... | |
Ifpack_IlukGraph (const Ifpack_IlukGraph &Graph_in) | |
Copy constructor. | |
virtual | ~Ifpack_IlukGraph () |
Ifpack_IlukGraph Destructor. | |
int | SetParameters (const Teuchos::ParameterList ¶meterlist, bool cerr_warning_if_unused=false) |
Set parameters using Teuchos::ParameterList object. | |
virtual int | ConstructFilledGraph () |
Does the actual construction of the graph. | |
virtual int | ConstructOverlapGraph () |
Does the actual construction of the overlap matrix graph. | |
virtual int | LevelFill () const |
Returns the level of fill used to construct this graph. | |
virtual int | LevelOverlap () const |
Returns the level of overlap used to construct this graph. | |
int | NumGlobalBlockRows () const |
Returns the number of global matrix rows. | |
int | NumGlobalBlockCols () const |
Returns the number of global matrix columns. | |
int | NumGlobalRows () const |
Returns the number of global matrix rows. | |
int | NumGlobalCols () const |
Returns the number of global matrix columns. | |
int | NumGlobalNonzeros () const |
Returns the number of nonzero entries in the global graph. | |
virtual int | NumGlobalBlockDiagonals () const |
Returns the number of diagonal entries found in the global input graph. | |
long long | NumGlobalBlockRows64 () const |
Returns the number of global matrix rows. | |
long long | NumGlobalBlockCols64 () const |
Returns the number of global matrix columns. | |
long long | NumGlobalRows64 () const |
Returns the number of global matrix rows. | |
long long | NumGlobalCols64 () const |
Returns the number of global matrix columns. | |
long long | NumGlobalNonzeros64 () const |
Returns the number of nonzero entries in the global graph. | |
virtual long long | NumGlobalBlockDiagonals64 () const |
Returns the number of diagonal entries found in the global input graph. | |
int | NumMyBlockRows () const |
Returns the number of local matrix rows. | |
int | NumMyBlockCols () const |
Returns the number of local matrix columns. | |
int | NumMyRows () const |
Returns the number of local matrix rows. | |
int | NumMyCols () const |
Returns the number of local matrix columns. | |
int | NumMyNonzeros () const |
Returns the number of nonzero entries in the local graph. | |
virtual int | NumMyBlockDiagonals () const |
Returns the number of diagonal entries found in the local input graph. | |
int | IndexBase () const |
Returns the index base for row and column indices for this graph. | |
long long | IndexBase64 () const |
virtual Epetra_CrsGraph & | L_Graph () |
Returns the graph of lower triangle of the ILU(k) graph as a Epetra_CrsGraph. | |
virtual Epetra_CrsGraph & | U_Graph () |
Returns the graph of upper triangle of the ILU(k) graph as a Epetra_CrsGraph. | |
virtual Epetra_CrsGraph & | L_Graph () const |
Returns the graph of lower triangle of the ILU(k) graph as a Epetra_CrsGraph. | |
virtual Epetra_CrsGraph & | U_Graph () const |
Returns the graph of upper triangle of the ILU(k) graph as a Epetra_CrsGraph. | |
virtual Epetra_Import * | OverlapImporter () const |
Returns the importer used to create the overlapped graph. | |
virtual Epetra_CrsGraph * | OverlapGraph () const |
Returns the the overlapped graph. | |
virtual const Epetra_BlockMap & | DomainMap () const |
Returns the Epetra_BlockMap object associated with the domain of this matrix operator. | |
virtual const Epetra_BlockMap & | RangeMap () const |
Returns the Epetra_BlockMap object associated with the range of this matrix operator. | |
virtual const Epetra_Comm & | Comm () const |
Returns the Epetra_BlockMap object associated with the range of this matrix operator. | |
Friends | |
std::ostream & | operator<< (std::ostream &os, const Ifpack_IlukGraph &A) |
<< operator will work for Ifpack_IlukGraph. | |
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class preconditioners.
The Ifpack_IlukGraph class enable the construction matrix graphs using level-fill algorithms. The only function required for construction is an ExtractRowView capability, i.e., the matrix that is passed in to the constructor must implement the Ifpack_CrsGraph interface defined in Ifpack_CrsMatrix.h
Constructing Ifpack_IlukGraph objects
Constructing Ifpack_IlukGraph objects is usually a two step process of passing in a Ifpack_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 Ifpack_IlukGraph object has two Epetra_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 Epetra_CrsGraph InsertIndices and RemoveIndices functions. However, in this case FillComplete must be called before the graph is used for subsequent operations.
Definition at line 77 of file Ifpack_IlukGraph.h.
Ifpack_IlukGraph::Ifpack_IlukGraph | ( | const Epetra_CrsGraph & | Graph_in, |
int | LevelFill_in, | ||
int | LevelOverlap_in | ||
) |
Ifpack_IlukGraph constuctor.
Creates a Ifpack_IlukGraph object using the input graph and specified level of fill.
In | Graph_in - An existing Ifpack_CrsGraph. This object must implement the Ifpack_CrsGraph functions that provide graph dimension and pattern information. |
In | LevelFill_in - The level of fill to compute via ILU(k) algorithm. |
In | LevelOverlap_in - The level of between subdomains. |
Definition at line 52 of file Ifpack_IlukGraph.cpp.