Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Ifpack_ILU.h
Go to the documentation of this file.
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_ILU_H
44 #define IFPACK_ILU_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 "Ifpack_Preconditioner.h"
54 #include "Ifpack_Condest.h"
55 #include "Ifpack_ScalingType.h"
56 #include "Ifpack_IlukGraph.h"
57 #include "Epetra_CompObject.h"
58 #include "Epetra_MultiVector.h"
59 #include "Epetra_Vector.h"
60 #include "Epetra_CrsGraph.h"
61 #include "Epetra_CrsMatrix.h"
62 #include "Epetra_BlockMap.h"
63 #include "Epetra_Map.h"
64 #include "Epetra_Object.h"
65 #include "Epetra_Comm.h"
66 #include "Epetra_RowMatrix.h"
67 #include "Epetra_Time.h"
68 #include "Teuchos_RefCountPtr.hpp"
69 
70 namespace Teuchos {
71  class ParameterList;
72 }
73 
75 
90 
91 public:
92  // @{ Constructors and destructors.
95 
98  {
99  Destroy();
100  }
101 
102  // @}
103  // @{ Construction methods
104 
106  int Initialize();
107 
109  bool IsInitialized() const
110  {
111  return(IsInitialized_);
112  }
113 
115 
123  int Compute();
124 
126  bool IsComputed() const
127  {
128  return(IsComputed_);
129  }
130 
132  /* This method is only available if the Teuchos package is enabled.
133  This method recognizes four parameter names: relax_value,
134  absolute_threshold, relative_threshold and overlap_mode. These names are
135  case insensitive, and in each case except overlap_mode, the ParameterEntry
136  must have type double. For overlap_mode, the ParameterEntry must have
137  type Epetra_CombineMode.
138  */
139  int SetParameters(Teuchos::ParameterList& parameterlist);
140 
142 
150  int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
151  // @}
152 
153  // @{ Mathematical functions.
154  // Applies the matrix to X, returns the result in Y.
155  int Apply(const Epetra_MultiVector& X,
156  Epetra_MultiVector& Y) const
157  {
158  return(Multiply(false,X,Y));
159  }
160 
161  int Multiply(bool Trans, const Epetra_MultiVector& X,
162  Epetra_MultiVector& Y) const;
163 
165 
178  int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
179 
181  double Condest(const Ifpack_CondestType CT = Ifpack_Cheap,
182  const int MaxIters = 1550,
183  const double Tol = 1e-9,
184  Epetra_RowMatrix* Matrix_in = 0);
185 
187  double Condest() const
188  {
189  return(Condest_);
190  }
191 
192  // @}
193  // @{ Query methods
194 
196  const Epetra_CrsMatrix & L() const {return(*L_);};
197 
199  const Epetra_Vector & D() const {return(*D_);};
200 
202  const Epetra_CrsMatrix & U() const {return(*U_);};
203 
205  const char* Label() const {return(Label_);}
206 
208  int SetLabel(const char* Label_in)
209  {
210  strcpy(Label_,Label_in);
211  return(0);
212  }
213 
215  double NormInf() const {return(0.0);};
216 
218  bool HasNormInf() const {return(false);};
219 
221  bool UseTranspose() const {return(UseTranspose_);};
222 
224  const Epetra_Map & OperatorDomainMap() const {return(U_->OperatorDomainMap());};
225 
227  const Epetra_Map & OperatorRangeMap() const{return(L_->OperatorRangeMap());};
228 
230  const Epetra_Comm & Comm() const{return(Comm_);};
231 
233  const Epetra_RowMatrix& Matrix() const
234  {
235  return(*A_);
236  }
237 
239  virtual std::ostream& Print(std::ostream& os) const;
240 
242  virtual int NumInitialize() const
243  {
244  return(NumInitialize_);
245  }
246 
248  virtual int NumCompute() const
249  {
250  return(NumCompute_);
251  }
252 
254  virtual int NumApplyInverse() const
255  {
256  return(NumApplyInverse_);
257  }
258 
260  virtual double InitializeTime() const
261  {
262  return(InitializeTime_);
263  }
264 
266  virtual double ComputeTime() const
267  {
268  return(ComputeTime_);
269  }
270 
272  virtual double ApplyInverseTime() const
273  {
274  return(ApplyInverseTime_);
275  }
276 
278  virtual double InitializeFlops() const
279  {
280  return(0.0);
281  }
282 
283  virtual double ComputeFlops() const
284  {
285  return(ComputeFlops_);
286  }
287 
288  virtual double ApplyInverseFlops() const
289  {
290  return(ApplyInverseFlops_);
291  }
292 
293 private:
294 
295  // @}
296  // @{ Private methods
297 
300  Comm_(RHS.Comm()),
301  Time_(RHS.Comm())
302  {}
303 
305  Ifpack_ILU& operator=(const Ifpack_ILU& /* RHS */)
306  {
307  return(*this);
308  }
309 
311  void Destroy();
312 
314 
324  int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
325 
326  int ComputeSetup();
327  int InitAllValues(const Epetra_RowMatrix & A, int MaxNumEntries);
328 
330  int LevelOfFill() const {return LevelOfFill_;}
331 
333  double RelaxValue() const {return RelaxValue_;}
334 
336  double AbsoluteThreshold() const {return Athresh_;}
337 
339  double RelativeThreshold() const {return Rthresh_;}
340 
341 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
342  int NumGlobalRows() const {return(Graph().NumGlobalRows());};
344 
346  int NumGlobalCols() const {return(Graph().NumGlobalCols());};
347 
349  int NumGlobalNonzeros() const {return(L().NumGlobalNonzeros()+U().NumGlobalNonzeros());};
350 
352  virtual int NumGlobalBlockDiagonals() const {return(Graph().NumGlobalBlockDiagonals());};
353 #endif
354 
355  long long NumGlobalRows64() const {return(Graph().NumGlobalRows64());};
356 
357  long long NumGlobalCols64() const {return(Graph().NumGlobalCols64());};
358 
359  long long NumGlobalNonzeros64() const {return(L().NumGlobalNonzeros64()+U().NumGlobalNonzeros64());};
360 
361  virtual long long NumGlobalBlockDiagonals64() const {return(Graph().NumGlobalBlockDiagonals64());};
362 
364  int NumMyRows() const {return(Graph().NumMyRows());};
365 
367  int NumMyCols() const {return(Graph().NumMyCols());};
368 
370  int NumMyNonzeros() const {return(L().NumMyNonzeros()+U().NumMyNonzeros());};
371 
373  virtual int NumMyBlockDiagonals() const {return(Graph().NumMyBlockDiagonals());};
374 
376  virtual int NumMyDiagonals() const {return(NumMyDiagonals_);};
377 
379 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
380  int IndexBase() const {return(Graph().IndexBase());};
381 #endif
382  long long IndexBase64() const {return(Graph().IndexBase64());};
383 
385  const Ifpack_IlukGraph & Graph() const {return(*Graph_);};
386 
389  {
390  return(*A_);
391  }
392 
393  // @}
394  // @{ Internal data
395 
397  Teuchos::RefCountPtr<Epetra_RowMatrix> A_;
398  Teuchos::RefCountPtr<Ifpack_IlukGraph> Graph_;
399  Teuchos::RefCountPtr<Epetra_CrsGraph> CrsGraph_;
400  Teuchos::RefCountPtr<Epetra_Map> IlukRowMap_;
401  Teuchos::RefCountPtr<Epetra_Map> IlukDomainMap_;
402  Teuchos::RefCountPtr<Epetra_Map> IlukRangeMap_;
407  Teuchos::RefCountPtr<Epetra_CrsMatrix> L_;
409  Teuchos::RefCountPtr<Epetra_CrsMatrix> U_;
410  Teuchos::RefCountPtr<Epetra_CrsGraph> L_Graph_;
411  Teuchos::RefCountPtr<Epetra_CrsGraph> U_Graph_;
413  Teuchos::RefCountPtr<Epetra_Vector> D_;
415 
419  bool Factored_;
421  double RelaxValue_;
423  double Athresh_;
425  double Rthresh_;
427  double Condest_;
435  char Label_[160];
441  mutable int NumApplyInverse_;
445  double ComputeTime_;
447  mutable double ApplyInverseTime_;
451  mutable double ApplyInverseFlops_;
454 
455 };
456 
457 #endif /* IFPACK_ILU_H */
char Label_[160]
Label of this object.
Definition: Ifpack_ILU.h:435
Teuchos::RefCountPtr< Epetra_CrsGraph > CrsGraph_
Definition: Ifpack_ILU.h:399
double InitializeTime_
Contains the time for all successful calls to Initialize().
Definition: Ifpack_ILU.h:443
Teuchos::RefCountPtr< Epetra_CrsGraph > L_Graph_
Definition: Ifpack_ILU.h:410
const Epetra_Comm & Comm_
Definition: Ifpack_ILU.h:405
const Ifpack_IlukGraph & Graph() const
Returns the address of the Ifpack_IlukGraph associated with this factored matrix. ...
Definition: Ifpack_ILU.h:385
double Athresh_
absolute threshold
Definition: Ifpack_ILU.h:423
virtual int NumMyDiagonals() const
Returns the number of nonzero diagonal values found in matrix.
Definition: Ifpack_ILU.h:376
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y...
Definition: Ifpack_ILU.cpp:575
virtual double InitializeFlops() const
Returns the number of flops in the initialization phase.
Definition: Ifpack_ILU.h:278
Teuchos::RefCountPtr< Epetra_Vector > D_
Diagonal of factors.
Definition: Ifpack_ILU.h:413
virtual int NumCompute() const
Returns the number of calls to Compute().
Definition: Ifpack_ILU.h:248
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
Definition: Ifpack_ILU.h:218
long long IndexBase64() const
Definition: Ifpack_ILU.h:382
virtual double ComputeFlops() const
Returns the number of flops in the computation phase.
Definition: Ifpack_ILU.h:283
int LevelOfFill() const
Returns the level of fill.
Definition: Ifpack_ILU.h:330
int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied.
Definition: Ifpack_ILU.h:150
bool ValuesInitialized_
Definition: Ifpack_ILU.h:418
long long NumGlobalCols64() const
Definition: Ifpack_ILU.h:357
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
Definition: Ifpack_ILU.h:224
int Multiply(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Definition: Ifpack_ILU.cpp:535
bool UseTranspose() const
Returns the current UseTranspose setting.
Definition: Ifpack_ILU.h:221
Epetra_RowMatrix & Matrix()
Returns a reference to the matrix.
Definition: Ifpack_ILU.h:388
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
Definition: Ifpack_ILU.h:349
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...
Ifpack_ILU(const Ifpack_ILU &RHS)
Copy constructor (should never be used)
Definition: Ifpack_ILU.h:299
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Definition: Ifpack_ILU.h:155
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
Definition: Ifpack_ILU.h:370
Teuchos::RefCountPtr< Ifpack_IlukGraph > Graph_
Definition: Ifpack_ILU.h:398
const Epetra_RowMatrix & Matrix() const
Returns a reference to the matrix to be preconditioned.
Definition: Ifpack_ILU.h:233
const Epetra_Map * U_DomainMap_
Definition: Ifpack_ILU.h:403
bool IsComputed_
If true, the preconditioner has been successfully computed.
Definition: Ifpack_ILU.h:433
Teuchos::RefCountPtr< Epetra_RowMatrix > A_
Pointer to the Epetra_RowMatrix to factorize.
Definition: Ifpack_ILU.h:397
const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
Definition: Ifpack_ILU.h:230
const Epetra_Map * L_RangeMap_
Definition: Ifpack_ILU.h:404
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
Definition: Ifpack_ILU.h:227
bool UseTranspose_
Definition: Ifpack_ILU.h:414
Ifpack_ILU: A class for constructing and using an incomplete lower/upper (ILU) factorization of a giv...
Definition: Ifpack_ILU.h:89
const Epetra_CrsMatrix & U() const
Returns the address of the L factor associated with this factored matrix.
Definition: Ifpack_ILU.h:202
virtual double ComputeTime() const
Returns the time spent in Compute().
Definition: Ifpack_ILU.h:266
Ifpack_CondestType
Ifpack_CondestType: enum to define the type of condition number estimate.
~Ifpack_ILU()
Destructor.
Definition: Ifpack_ILU.h:97
int SetLabel(const char *Label_in)
Sets label for this object.
Definition: Ifpack_ILU.h:208
bool IsComputed() const
If factor is completed, this query returns true, otherwise it returns false.
Definition: Ifpack_ILU.h:126
double ApplyInverseFlops_
Contain sthe number of flops for ApplyInverse().
Definition: Ifpack_ILU.h:451
int NumCompute_
Contains the number of successful call to Compute().
Definition: Ifpack_ILU.h:439
int Compute()
Compute ILU factors L and U using the specified graph, diagonal perturbation thresholds and relaxatio...
Definition: Ifpack_ILU.cpp:334
virtual double ApplyInverseTime() const
Returns the time spent in ApplyInverse().
Definition: Ifpack_ILU.h:272
Ifpack_ScalingType enumerable type.
int InitAllValues(const Epetra_RowMatrix &A, int MaxNumEntries)
int NumGlobalCols() const
Returns the number of global matrix columns.
Definition: Ifpack_ILU.h:346
int SetParameters(Teuchos::ParameterList &parameterlist)
Set parameters using a Teuchos::ParameterList object.
Definition: Ifpack_ILU.cpp:100
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.
int NumMyCols() const
Returns the number of local matrix columns.
Definition: Ifpack_ILU.h:367
Epetra_Time Time_
Used for timing issues.
Definition: Ifpack_ILU.h:453
long long NumGlobalRows64() const
Definition: Ifpack_ILU.h:355
Teuchos::RefCountPtr< Epetra_CrsMatrix > L_
Contains the L factors.
Definition: Ifpack_ILU.h:407
int NumGlobalRows() const
Returns the number of global matrix rows.
Definition: Ifpack_ILU.h:343
Teuchos::RefCountPtr< Epetra_Map > IlukRangeMap_
Definition: Ifpack_ILU.h:402
double AbsoluteThreshold() const
Get absolute threshold value.
Definition: Ifpack_ILU.h:336
int IndexBase() const
Returns the index base for row and column indices for this graph.
Definition: Ifpack_ILU.h:380
int Initialize()
Initialize the preconditioner, does not touch matrix values.
Definition: Ifpack_ILU.cpp:239
double ApplyInverseTime_
Contains the time for all successful calls to ApplyInverse().
Definition: Ifpack_ILU.h:447
const Epetra_CrsMatrix & L() const
Returns the address of the L factor associated with this factored matrix.
Definition: Ifpack_ILU.h:196
virtual long long NumGlobalBlockDiagonals64() const
Definition: Ifpack_ILU.h:361
bool Factored_
Definition: Ifpack_ILU.h:419
virtual double InitializeTime() const
Returns the time spent in Initialize().
Definition: Ifpack_ILU.h:260
Teuchos::RefCountPtr< Epetra_CrsMatrix > U_
Contains the U factors.
Definition: Ifpack_ILU.h:409
virtual int NumGlobalBlockDiagonals() const
Returns the number of diagonal entries found in the global input graph.
Definition: Ifpack_ILU.h:352
virtual double ApplyInverseFlops() const
Returns the number of flops in the application of the preconditioner.
Definition: Ifpack_ILU.h:288
int ComputeSetup()
Definition: Ifpack_ILU.cpp:115
bool IsInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack_ILU.h:109
cheap estimate
double Condest() const
Returns the computed estimated condition number, or -1.0 if not computed.
Definition: Ifpack_ILU.h:187
virtual std::ostream & Print(std::ostream &os) const
Prints on stream basic information about this object.
Definition: Ifpack_ILU.cpp:634
Ifpack_ILU(Epetra_RowMatrix *A)
Constructor.
Definition: Ifpack_ILU.cpp:65
int NumMyDiagonals_
Definition: Ifpack_ILU.h:416
int NumMyRows() const
Returns the number of local matrix rows.
Definition: Ifpack_ILU.h:364
long long NumGlobalNonzeros64() const
Definition: Ifpack_ILU.h:359
virtual int NumMyBlockDiagonals() const
Returns the number of diagonal entries found in the local input graph.
Definition: Ifpack_ILU.h:373
int NumApplyInverse_
Contains the number of successful call to ApplyInverse().
Definition: Ifpack_ILU.h:441
int NumInitialize_
Contains the number of successful calls to Initialize().
Definition: Ifpack_ILU.h:437
Teuchos::RefCountPtr< Epetra_Map > IlukRowMap_
Definition: Ifpack_ILU.h:400
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
Definition: Ifpack_ILU.h:215
void Destroy()
Destroys all internal data.
Definition: Ifpack_ILU.cpp:92
Teuchos::RefCountPtr< Epetra_CrsGraph > U_Graph_
Definition: Ifpack_ILU.h:411
Ifpack_ILU & operator=(const Ifpack_ILU &)
operator= (should never be used)
Definition: Ifpack_ILU.h:305
bool IsInitialized_
If true, the preconditioner has been successfully initialized.
Definition: Ifpack_ILU.h:431
bool Allocated_
Definition: Ifpack_ILU.h:417
Teuchos::RefCountPtr< Epetra_Map > IlukDomainMap_
Definition: Ifpack_ILU.h:401
double ComputeFlops_
Contains the number of flops for Compute().
Definition: Ifpack_ILU.h:449
const char * Label() const
Returns a character string describing the operator.
Definition: Ifpack_ILU.h:205
double RelativeThreshold() const
Get relative threshold value.
Definition: Ifpack_ILU.h:339
virtual int NumApplyInverse() const
Returns the number of calls to ApplyInverse().
Definition: Ifpack_ILU.h:254
virtual int NumInitialize() const
Returns the number of calls to Initialize().
Definition: Ifpack_ILU.h:242
double ComputeTime_
Contains the time for all successful calls to Compute().
Definition: Ifpack_ILU.h:445
const Epetra_Vector & D() const
Returns the address of the D factor associated with this factored matrix.
Definition: Ifpack_ILU.h:199
int Solve(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Ifpack_ILU forward/back solve on a Epetra_MultiVector X in Y...
Definition: Ifpack_ILU.cpp:500
double Rthresh_
relative threshold
Definition: Ifpack_ILU.h:425
double RelaxValue_
Relaxation value.
Definition: Ifpack_ILU.h:421
double RelaxValue() const
Get ILU(k) relaxation parameter.
Definition: Ifpack_ILU.h:333
#define RHS(a)
Definition: MatGenFD.c:60
int LevelOfFill_
Level of fill.
Definition: Ifpack_ILU.h:429
double Condest_
condition number estimate
Definition: Ifpack_ILU.h:427