43 #ifndef _IFPACK_ILUK_GRAPH_H_ 44 #define _IFPACK_ILUK_GRAPH_H_ 46 #include "Ifpack_ConfigDefs.h" 47 #include "Epetra_Object.h" 48 #include "Epetra_CrsGraph.h" 49 #include "Epetra_Import.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 140 return (
int)(NumGlobalBlockRows_);
142 throw "Ifpack_IlukGraph::NumGlobalBlockRows: GlobalIndices not int.";
148 return (
int)(NumGlobalBlockCols_);
150 throw "Ifpack_IlukGraph::NumGlobalBlockCols: GlobalIndices not int.";
156 return (
int)(NumGlobalRows_);
158 throw "Ifpack_IlukGraph::NumGlobalRows: GlobalIndices not int.";
164 return (
int)(NumGlobalCols_);
166 throw "Ifpack_IlukGraph::NumGlobalCols: GlobalIndices not int.";
171 return (
int)(NumGlobalNonzeros_);
173 throw "Ifpack_IlukGraph::NumGlobalNonzeros: GlobalIndices not int.";
179 return (
int)(NumGlobalBlockDiagonals_);
181 throw "Ifpack_IlukGraph::NumGlobalBlockDiagonals: GlobalIndices not int.";
222 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES 225 return (
int) IndexBase64();
226 throw "Ifpack_IlukGraph::IndexBase: GlobalIndices not int.";
229 long long IndexBase64()
const {
return(IndexBase_);};
265 Teuchos::RefCountPtr<Epetra_CrsGraph> OverlapGraph_;
266 Teuchos::RefCountPtr<Epetra_BlockMap> OverlapRowMap_;
267 Teuchos::RefCountPtr<Epetra_Import> OverlapImporter_;
270 Teuchos::RefCountPtr<Epetra_CrsGraph> L_Graph_;
271 Teuchos::RefCountPtr<Epetra_CrsGraph> U_Graph_;
272 long long IndexBase_;
273 long long NumGlobalRows_;
274 long long NumGlobalCols_;
275 long long NumGlobalBlockRows_;
276 long long NumGlobalBlockCols_;
277 long long NumGlobalBlockDiagonals_;
278 long long NumGlobalNonzeros_;
279 long long NumGlobalEntries_;
284 int NumMyBlockDiagonals_;
int NumGlobalBlockCols() const
Returns the number of global matrix columns.
int NumGlobalRows() const
Returns the number of global matrix rows.
virtual int NumGlobalBlockDiagonals() const
Returns the number of diagonal entries found in the global input graph.
int NumMyRows() const
Returns the number of local matrix rows.
int IndexBase() const
Returns the index base for row and column indices for this graph.
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.
virtual const Epetra_BlockMap & RangeMap() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
virtual Epetra_CrsGraph * OverlapGraph() const
Returns the the overlapped graph.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
virtual Epetra_CrsGraph & L_Graph()
Returns the graph of lower triangle of the ILU(k) graph as a Epetra_CrsGraph.
int SetParameters(const Teuchos::ParameterList ¶meterlist, bool cerr_warning_if_unused=false)
Set parameters using Teuchos::ParameterList object.
long long NumGlobalBlockCols64() const
Returns the number of global matrix columns.
virtual int LevelOverlap() const
Returns the level of overlap used to construct this graph.
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...
virtual int NumMyBlockDiagonals() const
Returns the number of diagonal entries found in the local input graph.
const Epetra_BlockMap & RowMap() const
int NumGlobalBlockRows() const
Returns the number of global matrix rows.
virtual int LevelFill() const
Returns the level of fill used to construct this graph.
virtual int ConstructOverlapGraph()
Does the actual construction of the overlap matrix graph.
virtual Epetra_CrsGraph & U_Graph()
Returns the graph of upper triangle of the ILU(k) graph as a Epetra_CrsGraph.
int NumGlobalCols() const
Returns the number of global matrix columns.
bool GlobalIndicesInt() const
long long NumGlobalRows64() const
Returns the number of global matrix rows.
virtual const Epetra_BlockMap & DomainMap() const
Returns the Epetra_BlockMap object associated with the domain of this matrix operator.
int NumMyBlockRows() const
Returns the number of local matrix rows.
long long NumGlobalBlockRows64() const
Returns the number of global matrix rows.
virtual const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
int NumMyBlockCols() const
Returns the number of local matrix columns.
virtual long long NumGlobalBlockDiagonals64() const
Returns the number of diagonal entries found in the global input graph.
long long NumGlobalCols64() const
Returns the number of global matrix columns.
virtual ~Ifpack_IlukGraph()
Ifpack_IlukGraph Destructor.
virtual Epetra_CrsGraph & U_Graph() const
Returns the graph of upper triangle of the ILU(k) graph as a Epetra_CrsGraph.
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
virtual int ConstructFilledGraph()
Does the actual construction of the graph.
int NumMyCols() const
Returns the number of local matrix columns.
virtual Epetra_CrsGraph & L_Graph() const
Returns the graph of lower triangle of the ILU(k) graph as a Epetra_CrsGraph.
friend std::ostream & operator<<(std::ostream &os, const Ifpack_IlukGraph &A)
<< operator will work for Ifpack_IlukGraph.
virtual Epetra_Import * OverlapImporter() const
Returns the importer used to create the overlapped graph.