IFPACK Development
Loading...
Searching...
No Matches
Ifpack_ILUT.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_ILUT_H
44#define IFPACK_ILUT_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_CondestType.h"
54#include "Ifpack_ScalingType.h"
55#include "Ifpack_Preconditioner.h"
56#include "Epetra_Vector.h"
57#include "Epetra_CrsMatrix.h"
58#include "Epetra_Time.h"
59#include "Teuchos_RefCountPtr.hpp"
60
63class Epetra_Comm;
64class Epetra_Map;
66
67namespace Teuchos {
68 class ParameterList;
69}
70
72
88
89public:
90 // @{ Constructors and Destructors
93
95 virtual ~Ifpack_ILUT();
96
97 // @}
98 // @{ Construction methods
100 /* This method is only available if the Teuchos package is enabled.
101 This method recognizes five parameter names: level_fill, drop_tolerance,
102 absolute_threshold, relative_threshold and overlap_mode. These names are
103 case insensitive. For level_fill the ParameterEntry must have type int, the
104 threshold entries must have type double and overlap_mode must have type
105 Epetra_CombineMode.
106 */
107 int SetParameters(Teuchos::ParameterList& parameterlis);
108
110
116 int Initialize();
117
119 bool IsInitialized() const
120 {
121 return(IsInitialized_);
122 }
123
125
133 int Compute();
134
136 bool IsComputed() const {return(IsComputed_);};
137
138 // Mathematical functions.
139
141
149 int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
150
151 int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
152
154 double Condest(const Ifpack_CondestType CT = Ifpack_Cheap,
155 const int MaxIters = 1550,
156 const double Tol = 1e-9,
157 Epetra_RowMatrix* Matrix_in = 0);
158
160 double Condest() const
161 {
162 return(Condest_);
163 }
164
166
174 int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
175
177 double NormInf() const {return(0.0);};
178
180 bool HasNormInf() const {return(false);};
181
183 bool UseTranspose() const {return(UseTranspose_);};
184
186 const Epetra_Map & OperatorDomainMap() const {return(A_.OperatorDomainMap());};
187
189 const Epetra_Map & OperatorRangeMap() const{return(A_.OperatorRangeMap());};
190
192 const Epetra_Comm & Comm() const{return(Comm_);};
193
196 {
197 return(A_);
198 }
199
201 const Epetra_CrsMatrix & L() const {return(*L_);};
202
204 const Epetra_CrsMatrix & U() const {return(*U_);};
205
207 const char* Label() const
208 {
209 return(Label_.c_str());
210 }
211
213 int SetLabel(const char* Label_in)
214 {
215 Label_ = Label_in;
216 return(0);
217 }
218
220 virtual std::ostream& Print(std::ostream& os) const;
221
223 virtual int NumInitialize() const
224 {
225 return(NumInitialize_);
226 }
227
229 virtual int NumCompute() const
230 {
231 return(NumCompute_);
232 }
233
235 virtual int NumApplyInverse() const
236 {
237 return(NumApplyInverse_);
238 }
239
241 virtual double InitializeTime() const
242 {
243 return(InitializeTime_);
244 }
245
247 virtual double ComputeTime() const
248 {
249 return(ComputeTime_);
250 }
251
253 virtual double ApplyInverseTime() const
254 {
255 return(ApplyInverseTime_);
256 }
257
259 virtual double InitializeFlops() const
260 {
261 return(0.0);
262 }
263
264 virtual double ComputeFlops() const
265 {
266 return(ComputeFlops_);
267 }
268
269 virtual double ApplyInverseFlops() const
270 {
271 return(ApplyInverseFlops_);
272 }
273
274 inline double LevelOfFill() const {
275 return(LevelOfFill_);
276 }
277
279 inline double RelaxValue() const {
280 return(Relax_);
281 }
282
284 inline double AbsoluteThreshold() const
285 {
286 return(Athresh_);
287 }
288
290 inline double RelativeThreshold() const
291 {
292 return(Rthresh_);
293 }
294
296 inline double DropTolerance() const
297 {
298 return(DropTolerance_);
299 }
300
302#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
303 int NumGlobalNonzeros() const {
304 // FIXME: diagonal of L_ should not be stored
305 return(L().NumGlobalNonzeros() + U().NumGlobalNonzeros() - L().NumGlobalRows());
306 }
307#endif
308 long long NumGlobalNonzeros64() const {
309 // FIXME: diagonal of L_ should not be stored
310 return(L().NumGlobalNonzeros64() + U().NumGlobalNonzeros64() - L().NumGlobalRows64());
311 }
312
314 int NumMyNonzeros() const {
315 return(L().NumMyNonzeros() + U().NumMyNonzeros());
316 }
317
318private:
319
320 // @}
321 // @{ Internal methods
322
324 Ifpack_ILUT(const Ifpack_ILUT& RHS) :
325 A_(RHS.Matrix()),
326 Comm_(RHS.Comm()),
327 Time_(Comm())
328 {};
329
331 Ifpack_ILUT& operator=(const Ifpack_ILUT& /* RHS */)
332 {
333 return(*this);
334 }
335
336 template<typename int_type>
337 int TCompute();
338
340 void Destroy();
341
342 // @}
343 // @{ Internal data
344
346 const Epetra_RowMatrix& A_;
348 const Epetra_Comm& Comm_;
350 Teuchos::RefCountPtr<Epetra_CrsMatrix> L_;
352 Teuchos::RefCountPtr<Epetra_CrsMatrix> U_;
354 double Condest_;
356 double Relax_;
358 double Athresh_;
360 double Rthresh_;
362 double LevelOfFill_;
364 double DropTolerance_;
366 std::string Label_;
368 bool IsInitialized_;
370 bool IsComputed_;
372 bool UseTranspose_;
374 int NumMyRows_;
376 int NumInitialize_;
378 int NumCompute_;
380 mutable int NumApplyInverse_;
382 double InitializeTime_;
384 double ComputeTime_;
386 mutable double ApplyInverseTime_;
388 double ComputeFlops_;
390 mutable double ApplyInverseFlops_;
392 mutable Epetra_Time Time_;
394 long long GlobalNonzeros_;
395 Teuchos::RefCountPtr<Epetra_SerialComm> SerialComm_;
396 Teuchos::RefCountPtr<Epetra_Map> SerialMap_;
397}; // Ifpack_ILUT
398
399#endif /* IFPACK_ILUT_H */
Ifpack_ScalingType enumerable type.
virtual const Epetra_Map & OperatorDomainMap() const=0
virtual const Epetra_Map & OperatorRangeMap() const=0
Ifpack_ILUT: A class for constructing and using an incomplete LU factorization of a given Epetra_RowM...
Definition Ifpack_ILUT.h:87
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
int Compute()
Compute IC factor U using the specified graph, diagonal perturbation thresholds and relaxation parame...
double RelativeThreshold() const
Get relative threshold value.
int SetParameters(Teuchos::ParameterList &parameterlis)
Set parameters using a Teuchos::ParameterList object.
double Condest() const
Returns the computed estimated condition number, or -1.0 if no computed.
virtual double InitializeTime() const
Returns the time spent in Initialize().
double DropTolerance() const
Gets the dropping tolerance.
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Ifpack_ILUT forward/back solve on a Epetra_MultiVector X in Y.
int SetLabel(const char *Label_in)
Sets the label for this object.
double RelaxValue() const
Set relative threshold value.
virtual int NumCompute() const
Returns the number of calls to Compute().
const Epetra_CrsMatrix & L() const
Returns a reference to the L factor.
const Epetra_CrsMatrix & U() const
Returns a reference to the U factor.
virtual int NumApplyInverse() const
Returns the number of calls to ApplyInverse().
virtual std::ostream & Print(std::ostream &os) const
Prints basic information on iostream. This function is used by operator<<.
bool IsInitialized() const
Returns true if the preconditioner has been successfully initialized.
const char * Label() const
Returns the label of this object.
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
virtual double ComputeFlops() const
Returns the number of flops in the computation phase.
virtual int NumInitialize() const
Returns the number of calls to Initialize().
int Initialize()
Initialize L and U with values from user matrix A.
virtual double ComputeTime() const
Returns the time spent in Compute().
const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
virtual ~Ifpack_ILUT()
Ifpack_ILUT Destructor.
int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied.
virtual double ApplyInverseTime() const
Returns the time spent in ApplyInverse().
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
bool IsComputed() const
If factor is completed, this query returns true, otherwise it returns false.
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
virtual double InitializeFlops() const
Returns the number of flops in the initialization phase.
bool UseTranspose() const
Returns the current UseTranspose setting.
virtual double ApplyInverseFlops() const
Returns the number of flops in the application of the preconditioner.
const Epetra_RowMatrix & Matrix() const
Returns a reference to the matrix to be preconditioned.
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
double AbsoluteThreshold() const
Get absolute threshold value.
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.