IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
Ifpack_IlukGraph.h
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack: Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2002) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #ifndef _IFPACK_ILUK_GRAPH_H_
44 #define _IFPACK_ILUK_GRAPH_H_
45 
46 #if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
47 #ifdef __GNUC__
48 #warning "The Ifpack package is deprecated"
49 #endif
50 #endif
51 
52 #include "Ifpack_ConfigDefs.h"
53 #include "Epetra_Object.h"
54 #include "Epetra_CrsGraph.h"
55 #include "Epetra_Import.h"
56 
57 #include "Teuchos_RefCountPtr.hpp"
58 
59 namespace Teuchos {
60  class ParameterList;
61 }
62 
64 
84 
85  // Give ostream << function some access to private and protected data/functions.
86 
87  friend std::ostream& operator << (std::ostream& os, const Ifpack_IlukGraph& A);
88 
89  public:
90 
92 
105  Ifpack_IlukGraph(const Epetra_CrsGraph & Graph_in, int LevelFill_in, int LevelOverlap_in);
106 
108  Ifpack_IlukGraph(const Ifpack_IlukGraph & Graph_in);
109 
111  virtual ~Ifpack_IlukGraph();
112 
114  /* This method is only available if the Teuchos package is enabled.
115  This method recogizes two parameter names: Level_fill and Level_overlap.
116  Both are case insensitive, and in both cases the ParameterEntry must
117  have type int.
118  */
119  int SetParameters(const Teuchos::ParameterList& parameterlist,
120  bool cerr_warning_if_unused=false);
121 
123  /*
124  \return Integer error code, set to 0 if successful.
125 
126  */
127  virtual int ConstructFilledGraph();
128 
130  /*
131  \return Integer error code, set to 0 if successful.
132 
133  */
134  virtual int ConstructOverlapGraph();
135 
137  virtual int LevelFill() const {return(LevelFill_);};
138 
140  virtual int LevelOverlap() const {return(LevelOverlap_);};
141 
142 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
143  int NumGlobalBlockRows() const {
145  if(Graph_.RowMap().GlobalIndicesInt())
146  return (int)(NumGlobalBlockRows_);
147  else
148  throw "Ifpack_IlukGraph::NumGlobalBlockRows: GlobalIndices not int.";
149  }
150 
152  int NumGlobalBlockCols() const {
153  if(Graph_.RowMap().GlobalIndicesInt())
154  return (int)(NumGlobalBlockCols_);
155  else
156  throw "Ifpack_IlukGraph::NumGlobalBlockCols: GlobalIndices not int.";
157  }
158 
160  int NumGlobalRows() const {
161  if(Graph_.RowMap().GlobalIndicesInt())
162  return (int)(NumGlobalRows_);
163  else
164  throw "Ifpack_IlukGraph::NumGlobalRows: GlobalIndices not int.";
165  }
166 
168  int NumGlobalCols() const {
169  if(Graph_.RowMap().GlobalIndicesInt())
170  return (int)(NumGlobalCols_);
171  else
172  throw "Ifpack_IlukGraph::NumGlobalCols: GlobalIndices not int.";
173  }
175  int NumGlobalNonzeros() const {
176  if(Graph_.RowMap().GlobalIndicesInt())
177  return (int)(NumGlobalNonzeros_);
178  else
179  throw "Ifpack_IlukGraph::NumGlobalNonzeros: GlobalIndices not int.";
180  }
181 
183  virtual int NumGlobalBlockDiagonals() const {
184  if(Graph_.RowMap().GlobalIndicesInt())
185  return (int)(NumGlobalBlockDiagonals_);
186  else
187  throw "Ifpack_IlukGraph::NumGlobalBlockDiagonals: GlobalIndices not int.";
188  }
189 #endif
190 
192  long long NumGlobalBlockRows64() const {return(NumGlobalBlockRows_);};
193 
195  long long NumGlobalBlockCols64() const {return(NumGlobalBlockCols_);};
196 
198  long long NumGlobalRows64() const {return(NumGlobalRows_);};
199 
201  long long NumGlobalCols64() const {return(NumGlobalCols_);};
203  long long NumGlobalNonzeros64() const {return(NumGlobalNonzeros_);};
204 
206  virtual long long NumGlobalBlockDiagonals64() const {return(NumGlobalBlockDiagonals_);};
207 
209  int NumMyBlockRows() const {return(NumMyBlockRows_);};
210 
212  int NumMyBlockCols() const {return(NumMyBlockCols_);};
213 
214 
216  int NumMyRows() const {return(NumMyRows_);};
217 
219  int NumMyCols() const {return(NumMyCols_);};
220 
222  int NumMyNonzeros() const {return(NumMyNonzeros_);};
223 
225  virtual int NumMyBlockDiagonals() const {return(NumMyBlockDiagonals_);};
226 
228 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
229  int IndexBase() const {
230  if(Graph_.RowMap().GlobalIndicesInt())
231  return (int) IndexBase64();
232  throw "Ifpack_IlukGraph::IndexBase: GlobalIndices not int.";
233  }
234 #endif
235  long long IndexBase64() const {return(IndexBase_);};
236 
238  virtual Epetra_CrsGraph & L_Graph() {return(*L_Graph_);};
239 
241  virtual Epetra_CrsGraph & U_Graph() {return(*U_Graph_);};
242 
244  virtual Epetra_CrsGraph & L_Graph() const {return(*L_Graph_);};
245 
247  virtual Epetra_CrsGraph & U_Graph() const {return(*U_Graph_);};
248 
250  virtual Epetra_Import * OverlapImporter() const {return(&*OverlapImporter_);};
251 
253  virtual Epetra_CrsGraph * OverlapGraph() const {return(&*OverlapGraph_);};
254 
256  virtual const Epetra_BlockMap & DomainMap() const {return(DomainMap_);};
257 
259  virtual const Epetra_BlockMap & RangeMap() const{return(RangeMap_);};
260 
262  virtual const Epetra_Comm & Comm() const{return(Comm_);};
263 
264  private:
265 
266 
267  const Epetra_CrsGraph & Graph_;
268  const Epetra_BlockMap & DomainMap_;
269  const Epetra_BlockMap & RangeMap_;
270  const Epetra_Comm & Comm_;
271  Teuchos::RefCountPtr<Epetra_CrsGraph> OverlapGraph_;
272  Teuchos::RefCountPtr<Epetra_BlockMap> OverlapRowMap_;
273  Teuchos::RefCountPtr<Epetra_Import> OverlapImporter_;
274  int LevelFill_;
275  int LevelOverlap_;
276  Teuchos::RefCountPtr<Epetra_CrsGraph> L_Graph_;
277  Teuchos::RefCountPtr<Epetra_CrsGraph> U_Graph_;
278  long long IndexBase_;
279  long long NumGlobalRows_;
280  long long NumGlobalCols_;
281  long long NumGlobalBlockRows_;
282  long long NumGlobalBlockCols_;
283  long long NumGlobalBlockDiagonals_;
284  long long NumGlobalNonzeros_;
285  long long NumGlobalEntries_;
286  int NumMyBlockRows_;
287  int NumMyBlockCols_;
288  int NumMyRows_;
289  int NumMyCols_;
290  int NumMyBlockDiagonals_;
291  int NumMyNonzeros_;
292  int NumMyEntries_;
293 
294 
295  };
296 
298 std::ostream& operator << (std::ostream& os, const Ifpack_IlukGraph& A);
299 
300 #endif /* _IFPACK_ILUK_GRAPH_H_ */
virtual const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
virtual Epetra_CrsGraph & L_Graph() const
Returns the graph of lower triangle of the ILU(k) graph as a Epetra_CrsGraph.
Ifpack_IlukGraph(const Epetra_CrsGraph &Graph_in, int LevelFill_in, int LevelOverlap_in)
Ifpack_IlukGraph constuctor.
int NumMyCols() 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.
int NumMyBlockCols() const
Returns the number of local matrix columns.
virtual Epetra_CrsGraph & L_Graph()
Returns the graph of lower triangle of the ILU(k) graph as a Epetra_CrsGraph.
virtual Epetra_CrsGraph * OverlapGraph() const
Returns the the overlapped graph.
int SetParameters(const Teuchos::ParameterList &parameterlist, bool cerr_warning_if_unused=false)
Set parameters using Teuchos::ParameterList object.
int NumMyRows() const
Returns the number of local matrix rows.
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...
bool GlobalIndicesInt() const
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
int NumGlobalBlockCols() const
Returns the number of global matrix columns.
virtual const Epetra_BlockMap & RangeMap() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
int IndexBase() const
Returns the index base for row and column indices for this graph.
virtual const Epetra_BlockMap & DomainMap() const
Returns the Epetra_BlockMap object associated with the domain of this matrix operator.
virtual int LevelFill() const
Returns the level of fill used to construct this graph.
long long NumGlobalBlockRows64() const
Returns the number of global matrix rows.
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.
virtual int NumGlobalBlockDiagonals() const
Returns the number of diagonal entries found in the global input graph.
long long NumGlobalCols64() const
Returns the number of global matrix columns.
long long NumGlobalRows64() const
Returns the number of global matrix rows.
const Epetra_BlockMap & RowMap() const
virtual int LevelOverlap() const
Returns the level of overlap used to construct this graph.
int NumGlobalRows() const
Returns the number of global matrix rows.
virtual Epetra_Import * OverlapImporter() const
Returns the importer used to create the overlapped graph.
long long NumGlobalBlockCols64() const
Returns the number of global matrix columns.
int NumGlobalCols() const
Returns the number of global matrix columns.
virtual int NumMyBlockDiagonals() const
Returns the number of diagonal entries found in the local input graph.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
virtual Epetra_CrsGraph & U_Graph() const
Returns the graph of upper triangle of the ILU(k) graph as a Epetra_CrsGraph.
int NumGlobalBlockRows() const
Returns the number of global matrix rows.
virtual ~Ifpack_IlukGraph()
Ifpack_IlukGraph Destructor.
virtual int ConstructFilledGraph()
Does the actual construction of the graph.
long long NumGlobalNonzeros64() const
Returns the number of nonzero entries in the global graph.
int NumMyBlockRows() const
Returns the number of local matrix rows.
friend std::ostream & operator<<(std::ostream &os, const Ifpack_IlukGraph &A)
&lt;&lt; operator will work for Ifpack_IlukGraph.