43#ifndef _IFPACK_ILUK_GRAPH_H_
44#define _IFPACK_ILUK_GRAPH_H_
51#include "Teuchos_RefCountPtr.hpp"
113 int SetParameters(
const Teuchos::ParameterList& parameterlist,
114 bool cerr_warning_if_unused=
false);
136#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
139 if(
Graph_.RowMap().GlobalIndicesInt())
142 throw "Ifpack_IlukGraph::NumGlobalBlockRows: GlobalIndices not int.";
147 if(
Graph_.RowMap().GlobalIndicesInt())
150 throw "Ifpack_IlukGraph::NumGlobalBlockCols: GlobalIndices not int.";
155 if(
Graph_.RowMap().GlobalIndicesInt())
158 throw "Ifpack_IlukGraph::NumGlobalRows: GlobalIndices not int.";
163 if(
Graph_.RowMap().GlobalIndicesInt())
166 throw "Ifpack_IlukGraph::NumGlobalCols: GlobalIndices not int.";
170 if(
Graph_.RowMap().GlobalIndicesInt())
173 throw "Ifpack_IlukGraph::NumGlobalNonzeros: GlobalIndices not int.";
178 if(
Graph_.RowMap().GlobalIndicesInt())
181 throw "Ifpack_IlukGraph::NumGlobalBlockDiagonals: GlobalIndices not int.";
222#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
224 if(
Graph_.RowMap().GlobalIndicesInt())
226 throw "Ifpack_IlukGraph::IndexBase: GlobalIndices not int.";
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 precondition...
virtual int LevelFill() const
Returns the level of fill used to construct this graph.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
int SetParameters(const Teuchos::ParameterList ¶meterlist, bool cerr_warning_if_unused=false)
Set parameters using Teuchos::ParameterList object.
long long NumGlobalNonzeros64() const
Returns the number of nonzero entries in the global graph.
Ifpack_IlukGraph(const Epetra_CrsGraph &Graph_in, int LevelFill_in, int LevelOverlap_in)
Ifpack_IlukGraph constuctor.
int NumGlobalBlockRows() const
Returns the number of global matrix rows.
virtual ~Ifpack_IlukGraph()
Ifpack_IlukGraph Destructor.
Teuchos::RefCountPtr< Epetra_BlockMap > OverlapRowMap_
virtual Epetra_CrsGraph & L_Graph() const
Returns the graph of lower triangle of the ILU(k) graph as a Epetra_CrsGraph.
int NumGlobalRows() const
Returns the number of global matrix rows.
long long NumGlobalNonzeros_
virtual Epetra_CrsGraph & U_Graph() const
Returns the graph of upper triangle of the ILU(k) graph as a Epetra_CrsGraph.
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.
const Epetra_BlockMap & DomainMap_
virtual const Epetra_BlockMap & DomainMap() const
Returns the Epetra_BlockMap object associated with the domain of this matrix operator.
virtual int ConstructOverlapGraph()
Does the actual construction of the overlap matrix graph.
long long NumGlobalEntries_
const Epetra_Comm & Comm_
const Epetra_CrsGraph & Graph_
virtual Epetra_CrsGraph & U_Graph()
Returns the graph of upper triangle of the ILU(k) graph as a Epetra_CrsGraph.
long long NumGlobalBlockRows64() const
Returns the number of global matrix rows.
int NumMyBlockRows() const
Returns the number of local matrix rows.
Teuchos::RefCountPtr< Epetra_Import > OverlapImporter_
Teuchos::RefCountPtr< Epetra_CrsGraph > OverlapGraph_
long long IndexBase64() const
int NumMyBlockCols() const
Returns the number of local matrix columns.
Teuchos::RefCountPtr< Epetra_CrsGraph > L_Graph_
int NumGlobalCols() const
Returns the number of global matrix columns.
int NumMyCols() const
Returns the number of local matrix columns.
long long NumGlobalBlockCols64() const
Returns the number of global matrix columns.
virtual long long NumGlobalBlockDiagonals64() const
Returns the number of diagonal entries found in the global input graph.
virtual int NumGlobalBlockDiagonals() const
Returns the number of diagonal entries found in the global input graph.
Teuchos::RefCountPtr< Epetra_CrsGraph > U_Graph_
virtual const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
const Epetra_BlockMap & RangeMap_
long long NumGlobalRows64() const
Returns the number of global matrix rows.
int NumMyRows() const
Returns the number of local matrix rows.
virtual Epetra_CrsGraph * OverlapGraph() const
Returns the the overlapped graph.
long long NumGlobalBlockRows_
virtual const Epetra_BlockMap & RangeMap() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
long long NumGlobalCols64() const
Returns the number of global matrix columns.
friend std::ostream & operator<<(std::ostream &os, const Ifpack_IlukGraph &A)
<< operator will work for Ifpack_IlukGraph.
virtual int LevelOverlap() const
Returns the level of overlap used to construct this graph.
virtual Epetra_CrsGraph & L_Graph()
Returns the graph of lower triangle of the ILU(k) graph as a Epetra_CrsGraph.
virtual int ConstructFilledGraph()
Does the actual construction of the graph.
long long NumGlobalBlockDiagonals_
long long NumGlobalBlockCols_
int NumGlobalBlockCols() const
Returns the number of global matrix columns.
virtual Epetra_Import * OverlapImporter() const
Returns the importer used to create the overlapped graph.