Ifpack Package Browser (Single Doxygen Collection)  Development
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 #include "Ifpack_ConfigDefs.h"
47 #include "Ifpack_Preconditioner.h"
48 #include "Ifpack_Condest.h"
49 #include "Ifpack_ScalingType.h"
50 #include "Ifpack_IlukGraph.h"
51 #include "Epetra_CompObject.h"
52 #include "Epetra_MultiVector.h"
53 #include "Epetra_Vector.h"
54 #include "Epetra_CrsGraph.h"
55 #include "Epetra_CrsMatrix.h"
56 #include "Epetra_BlockMap.h"
57 #include "Epetra_Map.h"
58 #include "Epetra_Object.h"
59 #include "Epetra_Comm.h"
60 #include "Epetra_RowMatrix.h"
61 #include "Epetra_Time.h"
62 #include "Teuchos_RefCountPtr.hpp"
63 
64 namespace Teuchos {
65  class ParameterList;
66 }
67 
69 
84 
85 public:
86  // @{ Constructors and destructors.
89 
92  {
93  Destroy();
94  }
95 
96  // @}
97  // @{ Construction methods
98 
100  int Initialize();
101 
103  bool IsInitialized() const
104  {
105  return(IsInitialized_);
106  }
107 
109 
117  int Compute();
118 
120  bool IsComputed() const
121  {
122  return(IsComputed_);
123  }
124 
126  /* This method is only available if the Teuchos package is enabled.
127  This method recognizes four parameter names: relax_value,
128  absolute_threshold, relative_threshold and overlap_mode. These names are
129  case insensitive, and in each case except overlap_mode, the ParameterEntry
130  must have type double. For overlap_mode, the ParameterEntry must have
131  type Epetra_CombineMode.
132  */
133  int SetParameters(Teuchos::ParameterList& parameterlist);
134 
136 
145  int SetUseTranspose(bool UseTranspose_in) {UseTranspose_ = UseTranspose_in; return(0);};
146  // @}
147 
148  // @{ Mathematical functions.
149  // Applies the matrix to X, returns the result in Y.
150  int Apply(const Epetra_MultiVector& X,
151  Epetra_MultiVector& Y) const
152  {
153  return(Multiply(false,X,Y));
154  }
155 
156  int Multiply(bool Trans, const Epetra_MultiVector& X,
157  Epetra_MultiVector& Y) const;
158 
160 
173  int ApplyInverse(const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
174 
176  double Condest(const Ifpack_CondestType CT = Ifpack_Cheap,
177  const int MaxIters = 1550,
178  const double Tol = 1e-9,
179  Epetra_RowMatrix* Matrix_in = 0);
180 
182  double Condest() const
183  {
184  return(Condest_);
185  }
186 
187  // @}
188  // @{ Query methods
189 
191  const Epetra_CrsMatrix & L() const {return(*L_);};
192 
194  const Epetra_Vector & D() const {return(*D_);};
195 
197  const Epetra_CrsMatrix & U() const {return(*U_);};
198 
200  const char* Label() const {return(Label_);}
201 
203  int SetLabel(const char* Label_in)
204  {
205  strcpy(Label_,Label_in);
206  return(0);
207  }
208 
210  double NormInf() const {return(0.0);};
211 
213  bool HasNormInf() const {return(false);};
214 
216  bool UseTranspose() const {return(UseTranspose_);};
217 
219  const Epetra_Map & OperatorDomainMap() const {return(U_->OperatorDomainMap());};
220 
222  const Epetra_Map & OperatorRangeMap() const{return(L_->OperatorRangeMap());};
223 
225  const Epetra_Comm & Comm() const{return(Comm_);};
226 
228  const Epetra_RowMatrix& Matrix() const
229  {
230  return(*A_);
231  }
232 
234  virtual std::ostream& Print(std::ostream& os) const;
235 
237  virtual int NumInitialize() const
238  {
239  return(NumInitialize_);
240  }
241 
243  virtual int NumCompute() const
244  {
245  return(NumCompute_);
246  }
247 
249  virtual int NumApplyInverse() const
250  {
251  return(NumApplyInverse_);
252  }
253 
255  virtual double InitializeTime() const
256  {
257  return(InitializeTime_);
258  }
259 
261  virtual double ComputeTime() const
262  {
263  return(ComputeTime_);
264  }
265 
267  virtual double ApplyInverseTime() const
268  {
269  return(ApplyInverseTime_);
270  }
271 
273  virtual double InitializeFlops() const
274  {
275  return(0.0);
276  }
277 
278  virtual double ComputeFlops() const
279  {
280  return(ComputeFlops_);
281  }
282 
283  virtual double ApplyInverseFlops() const
284  {
285  return(ApplyInverseFlops_);
286  }
287 
288 private:
289 
290  // @}
291  // @{ Private methods
292 
295  Comm_(RHS.Comm()),
296  Time_(RHS.Comm())
297  {}
298 
301  {
302  return(*this);
303  }
304 
306  void Destroy();
307 
309 
319  int Solve(bool Trans, const Epetra_MultiVector& X, Epetra_MultiVector& Y) const;
320 
321  int ComputeSetup();
322  int InitAllValues(const Epetra_RowMatrix & A, int MaxNumEntries);
323 
325  int LevelOfFill() const {return LevelOfFill_;}
326 
328  double RelaxValue() const {return RelaxValue_;}
329 
331  double AbsoluteThreshold() const {return Athresh_;}
332 
334  double RelativeThreshold() const {return Rthresh_;}
335 
336 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
337  int NumGlobalRows() const {return(Graph().NumGlobalRows());};
339 
341  int NumGlobalCols() const {return(Graph().NumGlobalCols());};
342 
344  int NumGlobalNonzeros() const {return(L().NumGlobalNonzeros()+U().NumGlobalNonzeros());};
345 
347  virtual int NumGlobalBlockDiagonals() const {return(Graph().NumGlobalBlockDiagonals());};
348 #endif
349 
350  long long NumGlobalRows64() const {return(Graph().NumGlobalRows64());};
351 
352  long long NumGlobalCols64() const {return(Graph().NumGlobalCols64());};
353 
354  long long NumGlobalNonzeros64() const {return(L().NumGlobalNonzeros64()+U().NumGlobalNonzeros64());};
355 
356  virtual long long NumGlobalBlockDiagonals64() const {return(Graph().NumGlobalBlockDiagonals64());};
357 
359  int NumMyRows() const {return(Graph().NumMyRows());};
360 
362  int NumMyCols() const {return(Graph().NumMyCols());};
363 
365  int NumMyNonzeros() const {return(L().NumMyNonzeros()+U().NumMyNonzeros());};
366 
368  virtual int NumMyBlockDiagonals() const {return(Graph().NumMyBlockDiagonals());};
369 
371  virtual int NumMyDiagonals() const {return(NumMyDiagonals_);};
372 
374 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
375  int IndexBase() const {return(Graph().IndexBase());};
376 #endif
377  long long IndexBase64() const {return(Graph().IndexBase64());};
378 
380  const Ifpack_IlukGraph & Graph() const {return(*Graph_);};
381 
384  {
385  return(*A_);
386  }
387 
388  // @}
389  // @{ Internal data
390 
392  Teuchos::RefCountPtr<Epetra_RowMatrix> A_;
393  Teuchos::RefCountPtr<Ifpack_IlukGraph> Graph_;
394  Teuchos::RefCountPtr<Epetra_CrsGraph> CrsGraph_;
395  Teuchos::RefCountPtr<Epetra_Map> IlukRowMap_;
396  Teuchos::RefCountPtr<Epetra_Map> IlukDomainMap_;
397  Teuchos::RefCountPtr<Epetra_Map> IlukRangeMap_;
402  Teuchos::RefCountPtr<Epetra_CrsMatrix> L_;
404  Teuchos::RefCountPtr<Epetra_CrsMatrix> U_;
405  Teuchos::RefCountPtr<Epetra_CrsGraph> L_Graph_;
406  Teuchos::RefCountPtr<Epetra_CrsGraph> U_Graph_;
408  Teuchos::RefCountPtr<Epetra_Vector> D_;
410 
414  bool Factored_;
416  double RelaxValue_;
418  double Athresh_;
420  double Rthresh_;
422  double Condest_;
430  char Label_[160];
436  mutable int NumApplyInverse_;
440  double ComputeTime_;
442  mutable double ApplyInverseTime_;
446  mutable double ApplyInverseFlops_;
449 
450 };
451 
452 #endif /* IFPACK_ILU_H */
Ifpack_ILU::OperatorDomainMap
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
Definition: Ifpack_ILU.h:219
Ifpack_ILU::Print
virtual std::ostream & Print(std::ostream &os) const
Prints on stream basic information about this object.
Definition: Ifpack_ILU.cpp:634
Ifpack_ILU::IlukDomainMap_
Teuchos::RefCountPtr< Epetra_Map > IlukDomainMap_
Definition: Ifpack_ILU.h:396
Ifpack_Preconditioner.h
Ifpack_IlukGraph.h
Ifpack_ILU::Ifpack_ILU
Ifpack_ILU(const Ifpack_ILU &RHS)
Copy constructor (should never be used)
Definition: Ifpack_ILU.h:294
Ifpack_ILU::InitializeTime
virtual double InitializeTime() const
Returns the time spent in Initialize().
Definition: Ifpack_ILU.h:255
Ifpack_ILU::IlukRowMap_
Teuchos::RefCountPtr< Epetra_Map > IlukRowMap_
Definition: Ifpack_ILU.h:395
Ifpack_ILU::U_DomainMap_
const Epetra_Map * U_DomainMap_
Definition: Ifpack_ILU.h:398
Ifpack_ILU::Destroy
void Destroy()
Destroys all internal data.
Definition: Ifpack_ILU.cpp:92
Ifpack_ILU::Graph
const Ifpack_IlukGraph & Graph() const
Returns the address of the Ifpack_IlukGraph associated with this factored matrix.
Definition: Ifpack_ILU.h:380
Ifpack_ILU::ComputeTime_
double ComputeTime_
Contains the time for all successful calls to Compute().
Definition: Ifpack_ILU.h:440
Ifpack_ILU::Matrix
Epetra_RowMatrix & Matrix()
Returns a reference to the matrix.
Definition: Ifpack_ILU.h:383
Ifpack_ILU::Condest
double Condest() const
Returns the computed estimated condition number, or -1.0 if not computed.
Definition: Ifpack_ILU.h:182
Ifpack_ILU::NumCompute
virtual int NumCompute() const
Returns the number of calls to Compute().
Definition: Ifpack_ILU.h:243
Ifpack_ILU::ComputeFlops_
double ComputeFlops_
Contains the number of flops for Compute().
Definition: Ifpack_ILU.h:444
Ifpack_ILU::UseTranspose_
bool UseTranspose_
Definition: Ifpack_ILU.h:409
Ifpack_ILU::Rthresh_
double Rthresh_
relative threshold
Definition: Ifpack_ILU.h:420
Ifpack_ILU::NumMyCols
int NumMyCols() const
Returns the number of local matrix columns.
Definition: Ifpack_ILU.h:362
RHS
#define RHS(a)
Definition: MatGenFD.c:60
Ifpack_ILU::NumGlobalRows
int NumGlobalRows() const
Returns the number of global matrix rows.
Definition: Ifpack_ILU.h:338
Ifpack_ILU::NumMyRows
int NumMyRows() const
Returns the number of local matrix rows.
Definition: Ifpack_ILU.h:359
Ifpack_ILU::OperatorRangeMap
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
Definition: Ifpack_ILU.h:222
Ifpack_ILU::SetLabel
int SetLabel(const char *Label_in)
Sets label for this object.
Definition: Ifpack_ILU.h:203
Ifpack_ILU::IsComputed_
bool IsComputed_
If true, the preconditioner has been successfully computed.
Definition: Ifpack_ILU.h:428
Ifpack_ILU::L_Graph_
Teuchos::RefCountPtr< Epetra_CrsGraph > L_Graph_
Definition: Ifpack_ILU.h:405
Ifpack_ILU::NumGlobalNonzeros64
long long NumGlobalNonzeros64() const
Definition: Ifpack_ILU.h:354
Ifpack_ILU::NumMyBlockDiagonals
virtual int NumMyBlockDiagonals() const
Returns the number of diagonal entries found in the local input graph.
Definition: Ifpack_ILU.h:368
Ifpack_ILU::Label
const char * Label() const
Returns a character string describing the operator.
Definition: Ifpack_ILU.h:200
Ifpack_ILU::Comm
const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
Definition: Ifpack_ILU.h:225
Ifpack_ILU::NumCompute_
int NumCompute_
Contains the number of successful call to Compute().
Definition: Ifpack_ILU.h:434
Ifpack_ILU::RelaxValue_
double RelaxValue_
Relaxation value.
Definition: Ifpack_ILU.h:416
Ifpack_ScalingType.h
Ifpack_ScalingType enumerable type.
Ifpack_ILU::Comm_
const Epetra_Comm & Comm_
Definition: Ifpack_ILU.h:400
Ifpack_ILU::IndexBase64
long long IndexBase64() const
Definition: Ifpack_ILU.h:377
Ifpack_ILU::SetUseTranspose
int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied.
Definition: Ifpack_ILU.h:145
Ifpack_ILU::IlukRangeMap_
Teuchos::RefCountPtr< Epetra_Map > IlukRangeMap_
Definition: Ifpack_ILU.h:397
Ifpack_IlukGraph
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...
Definition: Ifpack_IlukGraph.h:77
Ifpack_ILU::Athresh_
double Athresh_
absolute threshold
Definition: Ifpack_ILU.h:418
Ifpack_ILU::Initialize
int Initialize()
Initialize the preconditioner, does not touch matrix values.
Definition: Ifpack_ILU.cpp:239
Epetra_Comm
Ifpack_ILU::InitAllValues
int InitAllValues(const Epetra_RowMatrix &A, int MaxNumEntries)
Ifpack_ILU::IndexBase
int IndexBase() const
Returns the index base for row and column indices for this graph.
Definition: Ifpack_ILU.h:375
Ifpack_CondestType
Ifpack_CondestType
Ifpack_CondestType: enum to define the type of condition number estimate.
Definition: Ifpack_CondestType.h:48
Ifpack_ILU::L_
Teuchos::RefCountPtr< Epetra_CrsMatrix > L_
Contains the L factors.
Definition: Ifpack_ILU.h:402
Ifpack_Cheap
cheap estimate
Definition: Ifpack_CondestType.h:49
Ifpack_ILU::Factored_
bool Factored_
Definition: Ifpack_ILU.h:414
Ifpack_ILU::IsComputed
bool IsComputed() const
If factor is completed, this query returns true, otherwise it returns false.
Definition: Ifpack_ILU.h:120
Ifpack_ILU::RelaxValue
double RelaxValue() const
Get ILU(k) relaxation parameter.
Definition: Ifpack_ILU.h:328
Epetra_CrsMatrix
Ifpack_ILU::U_
Teuchos::RefCountPtr< Epetra_CrsMatrix > U_
Contains the U factors.
Definition: Ifpack_ILU.h:404
Ifpack_ILU::ComputeTime
virtual double ComputeTime() const
Returns the time spent in Compute().
Definition: Ifpack_ILU.h:261
Ifpack_ILU::ComputeSetup
int ComputeSetup()
Definition: Ifpack_ILU.cpp:115
Ifpack_ILU::L
const Epetra_CrsMatrix & L() const
Returns the address of the L factor associated with this factored matrix.
Definition: Ifpack_ILU.h:191
Ifpack_ILU::U
const Epetra_CrsMatrix & U() const
Returns the address of the L factor associated with this factored matrix.
Definition: Ifpack_ILU.h:197
Ifpack_ILU::Graph_
Teuchos::RefCountPtr< Ifpack_IlukGraph > Graph_
Definition: Ifpack_ILU.h:393
Ifpack_ILU::IsInitialized
bool IsInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack_ILU.h:103
Ifpack_ILU::Matrix
const Epetra_RowMatrix & Matrix() const
Returns a reference to the matrix to be preconditioned.
Definition: Ifpack_ILU.h:228
Ifpack_ConfigDefs.h
Ifpack_ILU::Ifpack_ILU
Ifpack_ILU(Epetra_RowMatrix *A)
Constructor.
Definition: Ifpack_ILU.cpp:65
Ifpack_ILU::operator=
Ifpack_ILU & operator=(const Ifpack_ILU &RHS)
operator= (should never be used)
Definition: Ifpack_ILU.h:300
Ifpack_ILU
Ifpack_ILU: A class for constructing and using an incomplete lower/upper (ILU) factorization of a giv...
Definition: Ifpack_ILU.h:83
Ifpack_ILU::NumGlobalCols64
long long NumGlobalCols64() const
Definition: Ifpack_ILU.h:352
Ifpack_ILU::NormInf
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
Definition: Ifpack_ILU.h:210
Epetra_RowMatrix
Ifpack_ILU::L_RangeMap_
const Epetra_Map * L_RangeMap_
Definition: Ifpack_ILU.h:399
Ifpack_ILU::RelativeThreshold
double RelativeThreshold() const
Get relative threshold value.
Definition: Ifpack_ILU.h:334
Ifpack_ILU::AbsoluteThreshold
double AbsoluteThreshold() const
Get absolute threshold value.
Definition: Ifpack_ILU.h:331
Ifpack_ILU::NumApplyInverse_
int NumApplyInverse_
Contains the number of successful call to ApplyInverse().
Definition: Ifpack_ILU.h:436
Ifpack_ILU::NumInitialize_
int NumInitialize_
Contains the number of successful calls to Initialize().
Definition: Ifpack_ILU.h:432
Ifpack_ILU::NumGlobalNonzeros
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
Definition: Ifpack_ILU.h:344
Ifpack_ILU::HasNormInf
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
Definition: Ifpack_ILU.h:213
Ifpack_ILU::Multiply
int Multiply(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Definition: Ifpack_ILU.cpp:535
Ifpack_ILU::A_
Teuchos::RefCountPtr< Epetra_RowMatrix > A_
Pointer to the Epetra_RowMatrix to factorize.
Definition: Ifpack_ILU.h:392
Ifpack_ILU::NumMyDiagonals
virtual int NumMyDiagonals() const
Returns the number of nonzero diagonal values found in matrix.
Definition: Ifpack_ILU.h:371
Ifpack_ILU::IsInitialized_
bool IsInitialized_
If true, the preconditioner has been successfully initialized.
Definition: Ifpack_ILU.h:426
Ifpack_ILU::NumGlobalRows64
long long NumGlobalRows64() const
Definition: Ifpack_ILU.h:350
Ifpack_ILU::NumApplyInverse
virtual int NumApplyInverse() const
Returns the number of calls to ApplyInverse().
Definition: Ifpack_ILU.h:249
Ifpack_ILU::Compute
int Compute()
Compute ILU factors L and U using the specified graph, diagonal perturbation thresholds and relaxatio...
Definition: Ifpack_ILU.cpp:334
Epetra_Vector
Ifpack_ILU::Apply
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Definition: Ifpack_ILU.h:150
Ifpack_ILU::SetParameters
int SetParameters(Teuchos::ParameterList &parameterlist)
Set parameters using a Teuchos::ParameterList object.
Definition: Ifpack_ILU.cpp:100
Ifpack_ILU::LevelOfFill
int LevelOfFill() const
Returns the level of fill.
Definition: Ifpack_ILU.h:325
Ifpack_ILU::U_Graph_
Teuchos::RefCountPtr< Epetra_CrsGraph > U_Graph_
Definition: Ifpack_ILU.h:406
Ifpack_ILU::~Ifpack_ILU
~Ifpack_ILU()
Destructor.
Definition: Ifpack_ILU.h:91
Ifpack_ILU::Allocated_
bool Allocated_
Definition: Ifpack_ILU.h:412
Epetra_Time
Ifpack_ILU::ApplyInverseFlops
virtual double ApplyInverseFlops() const
Returns the number of flops in the application of the preconditioner.
Definition: Ifpack_ILU.h:283
Ifpack_ILU::ApplyInverseFlops_
double ApplyInverseFlops_
Contain sthe number of flops for ApplyInverse().
Definition: Ifpack_ILU.h:446
Ifpack_ILU::D
const Epetra_Vector & D() const
Returns the address of the D factor associated with this factored matrix.
Definition: Ifpack_ILU.h:194
Ifpack_ILU::Condest_
double Condest_
condition number estimate
Definition: Ifpack_ILU.h:422
Ifpack_ILU::LevelOfFill_
int LevelOfFill_
Level of fill.
Definition: Ifpack_ILU.h:424
Ifpack_ILU::Label_
char Label_[160]
Label of this object.
Definition: Ifpack_ILU.h:430
Epetra_MultiVector
Ifpack_ILU::CrsGraph_
Teuchos::RefCountPtr< Epetra_CrsGraph > CrsGraph_
Definition: Ifpack_ILU.h:394
Ifpack_ILU::InitializeTime_
double InitializeTime_
Contains the time for all successful calls to Initialize().
Definition: Ifpack_ILU.h:438
A
Ifpack_ILU::Time_
Epetra_Time Time_
Used for timing issues.
Definition: Ifpack_ILU.h:448
Ifpack_Preconditioner
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.
Definition: Ifpack_Preconditioner.h:136
Ifpack_ILU::NumGlobalBlockDiagonals
virtual int NumGlobalBlockDiagonals() const
Returns the number of diagonal entries found in the global input graph.
Definition: Ifpack_ILU.h:347
Ifpack_ILU::NumGlobalCols
int NumGlobalCols() const
Returns the number of global matrix columns.
Definition: Ifpack_ILU.h:341
Ifpack_ILU::ApplyInverse
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
Teuchos::ParameterList
Ifpack_ILU::ApplyInverseTime_
double ApplyInverseTime_
Contains the time for all successful calls to ApplyInverse().
Definition: Ifpack_ILU.h:442
Teuchos
Ifpack_ILU::NumInitialize
virtual int NumInitialize() const
Returns the number of calls to Initialize().
Definition: Ifpack_ILU.h:237
Ifpack_ILU::ApplyInverseTime
virtual double ApplyInverseTime() const
Returns the time spent in ApplyInverse().
Definition: Ifpack_ILU.h:267
Ifpack_ILU::ComputeFlops
virtual double ComputeFlops() const
Returns the number of flops in the computation phase.
Definition: Ifpack_ILU.h:278
Ifpack_ILU::InitializeFlops
virtual double InitializeFlops() const
Returns the number of flops in the initialization phase.
Definition: Ifpack_ILU.h:273
Epetra_Map
Ifpack_ILU::D_
Teuchos::RefCountPtr< Epetra_Vector > D_
Diagonal of factors.
Definition: Ifpack_ILU.h:408
Ifpack_Condest.h
Ifpack_ILU::Solve
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
Ifpack_ILU::ValuesInitialized_
bool ValuesInitialized_
Definition: Ifpack_ILU.h:413
Ifpack_ILU::NumMyNonzeros
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
Definition: Ifpack_ILU.h:365
Ifpack_ILU::UseTranspose
bool UseTranspose() const
Returns the current UseTranspose setting.
Definition: Ifpack_ILU.h:216
Ifpack_ILU::NumGlobalBlockDiagonals64
virtual long long NumGlobalBlockDiagonals64() const
Definition: Ifpack_ILU.h:356
Ifpack_ILU::NumMyDiagonals_
int NumMyDiagonals_
Definition: Ifpack_ILU.h:411