IFPACK Development
Loading...
Searching...
No Matches
Ifpack_IKLU.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_IKLU_H
44#define IFPACK_IKLU_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 "Ifpack_IKLU_Utils.h"
57#include "Epetra_Vector.h"
58#include "Epetra_CrsMatrix.h"
59#include "Epetra_Time.h"
60#include "Teuchos_RefCountPtr.hpp"
61
64class Epetra_Comm;
65class Epetra_Map;
67
68namespace Teuchos {
69 class ParameterList;
70}
71
73
86
87public:
88 // @{ Constructors and Destructors
91
93 virtual ~Ifpack_IKLU();
94
95 // @}
96 // @{ Construction methods
98 /* This method is only available if the Teuchos package is enabled.
99 This method recognizes five parameter names: level_fill, drop_tolerance,
100 absolute_threshold, relative_threshold and overlap_mode. These names are
101 case insensitive. For level_fill the ParameterEntry must have type int, the
102 threshold entries must have type double and overlap_mode must have type
103 Epetra_CombineMode.
104 */
105 int SetParameters(Teuchos::ParameterList& parameterlis);
106
108
114 int Initialize();
115
117 bool IsInitialized() const
118 {
119 return(IsInitialized_);
120 }
121
123
131 int Compute();
132
134 bool IsComputed() const {return(IsComputed_);};
135
136 // Mathematical functions.
137
139
147 int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
148
149 int Apply(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
150
152 double Condest(const Ifpack_CondestType CT = Ifpack_Cheap,
153 const int MaxIters = 1550,
154 const double Tol = 1e-9,
155 Epetra_RowMatrix* Matrix_in = 0);
156
158 double Condest() const
159 {
160 return(Condest_);
161 }
162
164
173 int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
174
176 double NormInf() const {return(0.0);};
177
179 bool HasNormInf() const {return(false);};
180
182 bool UseTranspose() const {return(UseTranspose_);};
183
185 const Epetra_Map & OperatorDomainMap() const {return(A_.OperatorDomainMap());};
186
188 const Epetra_Map & OperatorRangeMap() const{return(A_.OperatorRangeMap());};
189
191 const Epetra_Comm & Comm() const{return(Comm_);};
192
195 {
196 return(A_);
197 }
198
200 const Epetra_CrsMatrix & L() const {return(*L_);};
201
203 const Epetra_CrsMatrix & U() const {return(*U_);};
204
206 const char* Label() const
207 {
208 return(Label_.c_str());
209 }
210
212 int SetLabel(const char* Label_in)
213 {
214 Label_ = Label_in;
215 return(0);
216 }
217
219 virtual std::ostream& Print(std::ostream& os) const;
220
222 virtual int NumInitialize() const
223 {
224 return(NumInitialize_);
225 }
226
228 virtual int NumCompute() const
229 {
230 return(NumCompute_);
231 }
232
234 virtual int NumApplyInverse() const
235 {
236 return(NumApplyInverse_);
237 }
238
240 virtual double InitializeTime() const
241 {
242 return(InitializeTime_);
243 }
244
246 virtual double ComputeTime() const
247 {
248 return(ComputeTime_);
249 }
250
252 virtual double ApplyInverseTime() const
253 {
254 return(ApplyInverseTime_);
255 }
256
258 virtual double InitializeFlops() const
259 {
260 return(0.0);
261 }
262
263 virtual double ComputeFlops() const
264 {
265 return(ComputeFlops_);
266 }
267
268 virtual double ApplyInverseFlops() const
269 {
270 return(ApplyInverseFlops_);
271 }
272
273 inline double LevelOfFill() const {
274 return(LevelOfFill_);
275 }
276
278 inline double RelaxValue() const {
279 return(Relax_);
280 }
281
283 inline double AbsoluteThreshold() const
284 {
285 return(Athresh_);
286 }
287
289 inline double RelativeThreshold() const
290 {
291 return(Rthresh_);
292 }
293
295 inline double DropTolerance() const
296 {
297 return(DropTolerance_);
298 }
299
301#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
302 int NumGlobalNonzeros() const {
303 // FIXME: diagonal of L_ should not be stored
304 return(L().NumGlobalNonzeros() + U().NumGlobalNonzeros() - L().NumGlobalRows());
305 }
306#endif
307 long long NumGlobalNonzeros64() const {
308 // FIXME: diagonal of L_ should not be stored
309 return(L().NumGlobalNonzeros64() + U().NumGlobalNonzeros64() - L().NumGlobalRows64());
310 }
311
313 int NumMyNonzeros() const {
314 return(L().NumMyNonzeros() + U().NumMyNonzeros());
315 }
316
317private:
318
319 // @}
320 // @{ Internal methods
321
323 Ifpack_IKLU(const Ifpack_IKLU& RHS) :
324 A_(RHS.Matrix()),
325 Comm_(RHS.Comm()),
326 Time_(Comm())
327 {};
328
330 Ifpack_IKLU& operator=(const Ifpack_IKLU& /* RHS */)
331 {
332 return(*this);
333 }
334
336 void Destroy();
337
338 // @}
339 // @{ Internal data
340
342 const Epetra_RowMatrix& A_;
344 const Epetra_Comm& Comm_;
346 Teuchos::RefCountPtr<Epetra_CrsMatrix> L_;
348 Teuchos::RefCountPtr<Epetra_CrsMatrix> U_;
350 double Condest_;
352 double Relax_;
354 double Athresh_;
356 double Rthresh_;
358 double LevelOfFill_;
360 double DropTolerance_;
362 std::string Label_;
364 bool IsInitialized_;
366 bool IsComputed_;
368 bool UseTranspose_;
370 int NumMyRows_;
372 int NumMyNonzeros_;
374 int NumInitialize_;
376 int NumCompute_;
378 mutable int NumApplyInverse_;
380 double InitializeTime_;
382 double ComputeTime_;
384 mutable double ApplyInverseTime_;
386 double ComputeFlops_;
388 mutable double ApplyInverseFlops_;
390 mutable Epetra_Time Time_;
392 long long GlobalNonzeros_;
393 Teuchos::RefCountPtr<Epetra_SerialComm> SerialComm_;
394 Teuchos::RefCountPtr<Epetra_Map> SerialMap_;
395
397 csr* csrA_;
398 css* cssS_;
400 csrn* csrnN_;
401
402}; // Ifpack_IKLU
403
404#endif /* IFPACK_IKLU_H */
Ifpack_ScalingType enumerable type.
virtual const Epetra_Map & OperatorDomainMap() const=0
virtual const Epetra_Map & OperatorRangeMap() const=0
Ifpack_IKLU: A class for constructing and using an incomplete LU factorization of a given Epetra_RowM...
Definition Ifpack_IKLU.h:85
virtual double ApplyInverseTime() const
Returns the time spent in ApplyInverse().
bool IsInitialized() const
Returns true if the preconditioner has been successfully initialized.
virtual int NumCompute() const
Returns the number of calls to Compute().
double DropTolerance() const
Gets the dropping tolerance.
virtual ~Ifpack_IKLU()
Ifpack_IKLU Destructor.
int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied.
bool IsComputed() const
If factor is completed, this query returns true, otherwise it returns false.
int Compute()
Compute IC factor U using the specified graph, diagonal perturbation thresholds and relaxation parame...
int SetParameters(Teuchos::ParameterList &parameterlis)
Set parameters using a Teuchos::ParameterList object.
const Epetra_RowMatrix & Matrix() const
Returns a reference to the matrix to be preconditioned.
const char * Label() const
Returns the label of this object.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
virtual int NumInitialize() const
Returns the number of calls to Initialize().
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
double AbsoluteThreshold() const
Get absolute threshold value.
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
virtual double ComputeFlops() const
Returns the number of flops in the computation phase.
double RelativeThreshold() const
Get relative threshold value.
virtual double InitializeTime() const
Returns the time spent in Initialize().
bool UseTranspose() const
Returns the current UseTranspose setting.
double Condest() const
Returns the computed estimated condition number, or -1.0 if no computed.
virtual double ComputeTime() const
Returns the time spent in Compute().
int Initialize()
Initialize L and U with values from user matrix A.
const Epetra_CrsMatrix & L() const
Returns a reference to the L factor.
const Epetra_CrsMatrix & U() const
Returns a reference to the U factor.
double RelaxValue() const
Set relative threshold value.
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
virtual int NumApplyInverse() const
Returns the number of calls to ApplyInverse().
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Ifpack_IKLU forward/back solve on a Epetra_MultiVector X in Y.
virtual double InitializeFlops() const
Returns the number of flops in the initialization phase.
const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
virtual std::ostream & Print(std::ostream &os) const
Prints basic information on iostream. This function is used by operator<<.
int SetLabel(const char *Label_in)
Sets the label for this object.
virtual double ApplyInverseFlops() const
Returns the number of flops in the application of the preconditioner.
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.