Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
Ifpack_ILU Class Reference

Ifpack_ILU: A class for constructing and using an incomplete lower/upper (ILU) factorization of a given Epetra_RowMatrix. More...

#include <Ifpack_ILU.h>

Inheritance diagram for Ifpack_ILU:
Inheritance graph
[legend]
 Ifpack_ILU (Epetra_RowMatrix *A)
 Constructor. More...
 
 ~Ifpack_ILU ()
 Destructor. More...
 
int Initialize ()
 Initialize the preconditioner, does not touch matrix values. More...
 
bool IsInitialized () const
 Returns true if the preconditioner has been successfully initialized. More...
 
int Compute ()
 Compute ILU factors L and U using the specified graph, diagonal perturbation thresholds and relaxation parameters. More...
 
bool IsComputed () const
 If factor is completed, this query returns true, otherwise it returns false. More...
 
int SetParameters (Teuchos::ParameterList &parameterlist)
 Set parameters using a Teuchos::ParameterList object. More...
 
int SetUseTranspose (bool UseTranspose_in)
 If set true, transpose of this operator will be applied. More...
 
int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 
int Multiply (bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 
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. More...
 
double Condest (const Ifpack_CondestType CT=Ifpack_Cheap, const int MaxIters=1550, const double Tol=1e-9, Epetra_RowMatrix *Matrix_in=0)
 Computes the estimated condition number and returns the value. More...
 
double Condest () const
 Returns the computed estimated condition number, or -1.0 if not computed. More...
 
const Epetra_CrsMatrixL () const
 Returns the address of the L factor associated with this factored matrix. More...
 
const Epetra_VectorD () const
 Returns the address of the D factor associated with this factored matrix. More...
 
const Epetra_CrsMatrixU () const
 Returns the address of the L factor associated with this factored matrix. More...
 
const char * Label () const
 Returns a character string describing the operator. More...
 
int SetLabel (const char *Label_in)
 Sets label for this object. More...
 
double NormInf () const
 Returns 0.0 because this class cannot compute Inf-norm. More...
 
bool HasNormInf () const
 Returns false because this class cannot compute an Inf-norm. More...
 
bool UseTranspose () const
 Returns the current UseTranspose setting. More...
 
const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this operator. More...
 
const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this operator. More...
 
const Epetra_CommComm () const
 Returns the Epetra_BlockMap object associated with the range of this matrix operator. More...
 
const Epetra_RowMatrixMatrix () const
 Returns a reference to the matrix to be preconditioned. More...
 
virtual std::ostream & Print (std::ostream &os) const
 Prints on stream basic information about this object. More...
 
virtual int NumInitialize () const
 Returns the number of calls to Initialize(). More...
 
virtual int NumCompute () const
 Returns the number of calls to Compute(). More...
 
virtual int NumApplyInverse () const
 Returns the number of calls to ApplyInverse(). More...
 
virtual double InitializeTime () const
 Returns the time spent in Initialize(). More...
 
virtual double ComputeTime () const
 Returns the time spent in Compute(). More...
 
virtual double ApplyInverseTime () const
 Returns the time spent in ApplyInverse(). More...
 
virtual double InitializeFlops () const
 Returns the number of flops in the initialization phase. More...
 
virtual double ComputeFlops () const
 Returns the number of flops in the computation phase. More...
 
virtual double ApplyInverseFlops () const
 Returns the number of flops in the application of the preconditioner. More...
 
 Ifpack_ILU (const Ifpack_ILU &RHS)
 Copy constructor (should never be used) More...
 
Ifpack_ILUoperator= (const Ifpack_ILU &)
 operator= (should never be used) More...
 
void Destroy ()
 Destroys all internal data. More...
 
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. More...
 
int ComputeSetup ()
 
int InitAllValues (const Epetra_RowMatrix &A, int MaxNumEntries)
 
int LevelOfFill () const
 Returns the level of fill. More...
 
double RelaxValue () const
 Get ILU(k) relaxation parameter. More...
 
double AbsoluteThreshold () const
 Get absolute threshold value. More...
 
double RelativeThreshold () const
 Get relative threshold value. More...
 
int NumGlobalRows () const
 Returns the number of global matrix rows. More...
 
int NumGlobalCols () const
 Returns the number of global matrix columns. More...
 
int NumGlobalNonzeros () const
 Returns the number of nonzero entries in the global graph. More...
 
virtual int NumGlobalBlockDiagonals () const
 Returns the number of diagonal entries found in the global input graph. More...
 
long long NumGlobalRows64 () const
 
long long NumGlobalCols64 () const
 
long long NumGlobalNonzeros64 () const
 
virtual long long NumGlobalBlockDiagonals64 () const
 
int NumMyRows () const
 Returns the number of local matrix rows. More...
 
int NumMyCols () const
 Returns the number of local matrix columns. More...
 
int NumMyNonzeros () const
 Returns the number of nonzero entries in the local graph. More...
 
virtual int NumMyBlockDiagonals () const
 Returns the number of diagonal entries found in the local input graph. More...
 
virtual int NumMyDiagonals () const
 Returns the number of nonzero diagonal values found in matrix. More...
 
int IndexBase () const
 Returns the index base for row and column indices for this graph. More...
 
long long IndexBase64 () const
 
const Ifpack_IlukGraphGraph () const
 Returns the address of the Ifpack_IlukGraph associated with this factored matrix. More...
 
Epetra_RowMatrixMatrix ()
 Returns a reference to the matrix. More...
 
Teuchos::RefCountPtr
< Epetra_RowMatrix
A_
 Pointer to the Epetra_RowMatrix to factorize. More...
 
Teuchos::RefCountPtr
< Ifpack_IlukGraph
Graph_
 
Teuchos::RefCountPtr
< Epetra_CrsGraph
CrsGraph_
 
Teuchos::RefCountPtr< Epetra_MapIlukRowMap_
 
Teuchos::RefCountPtr< Epetra_MapIlukDomainMap_
 
Teuchos::RefCountPtr< Epetra_MapIlukRangeMap_
 
const Epetra_MapU_DomainMap_
 
const Epetra_MapL_RangeMap_
 
const Epetra_CommComm_
 
Teuchos::RefCountPtr
< Epetra_CrsMatrix
L_
 Contains the L factors. More...
 
Teuchos::RefCountPtr
< Epetra_CrsMatrix
U_
 Contains the U factors. More...
 
Teuchos::RefCountPtr
< Epetra_CrsGraph
L_Graph_
 
Teuchos::RefCountPtr
< Epetra_CrsGraph
U_Graph_
 
Teuchos::RefCountPtr
< Epetra_Vector
D_
 Diagonal of factors. More...
 
bool UseTranspose_
 
int NumMyDiagonals_
 
bool Allocated_
 
bool ValuesInitialized_
 
bool Factored_
 
double RelaxValue_
 Relaxation value. More...
 
double Athresh_
 absolute threshold More...
 
double Rthresh_
 relative threshold More...
 
double Condest_
 condition number estimate More...
 
int LevelOfFill_
 Level of fill. More...
 
bool IsInitialized_
 If true, the preconditioner has been successfully initialized. More...
 
bool IsComputed_
 If true, the preconditioner has been successfully computed. More...
 
char Label_ [160]
 Label of this object. More...
 
int NumInitialize_
 Contains the number of successful calls to Initialize(). More...
 
int NumCompute_
 Contains the number of successful call to Compute(). More...
 
int NumApplyInverse_
 Contains the number of successful call to ApplyInverse(). More...
 
double InitializeTime_
 Contains the time for all successful calls to Initialize(). More...
 
double ComputeTime_
 Contains the time for all successful calls to Compute(). More...
 
double ApplyInverseTime_
 Contains the time for all successful calls to ApplyInverse(). More...
 
double ComputeFlops_
 Contains the number of flops for Compute(). More...
 
double ApplyInverseFlops_
 Contain sthe number of flops for ApplyInverse(). More...
 
Epetra_Time Time_
 Used for timing issues. More...
 

Additional Inherited Members

Detailed Description

Ifpack_ILU: A class for constructing and using an incomplete lower/upper (ILU) factorization of a given Epetra_RowMatrix.

The Ifpack_ILU class computes a "relaxed" ILU factorization with level k fill of a given Epetra_RowMatrix. The notion of relaxation is the same as described in "Experimental study of ILU preconditioners for indefinite matrices" by Chow and Saad.

Please refer to ifp_ilu for a general description of the ILU algorithm.

The complete list of supported parameters is reported in page ifp_params.

Author
Mike Heroux, Marzio Sala, SNL 9214.
Date
Last modified on 22-Jan-05.

Definition at line 83 of file Ifpack_ILU.h.

Constructor & Destructor Documentation

Ifpack_ILU::Ifpack_ILU ( Epetra_RowMatrix A)

Constructor.

Definition at line 65 of file Ifpack_ILU.cpp.

Ifpack_ILU::~Ifpack_ILU ( )
inline

Destructor.

Definition at line 91 of file Ifpack_ILU.h.

Ifpack_ILU::Ifpack_ILU ( const Ifpack_ILU RHS)
inlineprivate

Copy constructor (should never be used)

Definition at line 293 of file Ifpack_ILU.h.

Member Function Documentation

int Ifpack_ILU::Initialize ( )
virtual

Initialize the preconditioner, does not touch matrix values.

Implements Ifpack_Preconditioner.

Definition at line 239 of file Ifpack_ILU.cpp.

bool Ifpack_ILU::IsInitialized ( ) const
inlinevirtual

Returns true if the preconditioner has been successfully initialized.

Implements Ifpack_Preconditioner.

Definition at line 103 of file Ifpack_ILU.h.

int Ifpack_ILU::Compute ( )
virtual

Compute ILU factors L and U using the specified graph, diagonal perturbation thresholds and relaxation parameters.

This function computes the ILU(k) factors L and U using the current:

  1. Ifpack_IlukGraph specifying the structure of L and U.
  2. Value for the ILU(k) relaxation parameter.
  3. Value for the a priori diagonal threshold values.

InitValues() must be called before the factorization can proceed.

Implements Ifpack_Preconditioner.

Definition at line 334 of file Ifpack_ILU.cpp.

bool Ifpack_ILU::IsComputed ( ) const
inlinevirtual

If factor is completed, this query returns true, otherwise it returns false.

Implements Ifpack_Preconditioner.

Definition at line 120 of file Ifpack_ILU.h.

int Ifpack_ILU::SetParameters ( Teuchos::ParameterList parameterlist)
virtual

Set parameters using a Teuchos::ParameterList object.

Implements Ifpack_Preconditioner.

Definition at line 100 of file Ifpack_ILU.cpp.

int Ifpack_ILU::SetUseTranspose ( bool  UseTranspose_in)
inlinevirtual

If set true, transpose of this operator will be applied.

This flag allows the transpose of the given operator to be used implicitly. Setting this flag affects only the Apply() and ApplyInverse() methods. If the implementation of this interface does not support transpose use, this method should return a value of -1.

Parameters
UseTranspose_in- (In) If true, multiply by the transpose of operator, otherwise just use operator.
Returns
Always returns 0.

Implements Epetra_Operator.

Definition at line 144 of file Ifpack_ILU.h.

int Ifpack_ILU::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

Implements Epetra_Operator.

Definition at line 149 of file Ifpack_ILU.h.

int Ifpack_ILU::Multiply ( bool  Trans,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const

Definition at line 535 of file Ifpack_ILU.cpp.

int Ifpack_ILU::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.

In this implementation, we use several existing attributes to determine how virtual method ApplyInverse() should call the concrete method Solve(). We pass in the UpperTriangular(), the Epetra_CrsMatrix::UseTranspose(), and NoDiagonal() methods. The most notable warning is that if a matrix has no diagonal values we assume that there is an implicit unit diagonal that should be accounted for when doing a triangular solve.

Parameters
X- (In) A Epetra_MultiVector of dimension NumVectors to solve for.
OutY - (Out) A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Implements Ifpack_Preconditioner.

Definition at line 575 of file Ifpack_ILU.cpp.

double Ifpack_ILU::Condest ( const Ifpack_CondestType  CT = Ifpack_Cheap,
const int  MaxIters = 1550,
const double  Tol = 1e-9,
Epetra_RowMatrix Matrix_in = 0 
)
virtual

Computes the estimated condition number and returns the value.

Implements Ifpack_Preconditioner.

Definition at line 615 of file Ifpack_ILU.cpp.

double Ifpack_ILU::Condest ( ) const
inlinevirtual

Returns the computed estimated condition number, or -1.0 if not computed.

Implements Ifpack_Preconditioner.

Definition at line 181 of file Ifpack_ILU.h.

const Epetra_CrsMatrix& Ifpack_ILU::L ( ) const
inline

Returns the address of the L factor associated with this factored matrix.

Definition at line 190 of file Ifpack_ILU.h.

const Epetra_Vector& Ifpack_ILU::D ( ) const
inline

Returns the address of the D factor associated with this factored matrix.

Definition at line 193 of file Ifpack_ILU.h.

const Epetra_CrsMatrix& Ifpack_ILU::U ( ) const
inline

Returns the address of the L factor associated with this factored matrix.

Definition at line 196 of file Ifpack_ILU.h.

const char* Ifpack_ILU::Label ( ) const
inlinevirtual

Returns a character string describing the operator.

Implements Epetra_Operator.

Definition at line 199 of file Ifpack_ILU.h.

int Ifpack_ILU::SetLabel ( const char *  Label_in)
inline

Sets label for this object.

Definition at line 202 of file Ifpack_ILU.h.

double Ifpack_ILU::NormInf ( ) const
inlinevirtual

Returns 0.0 because this class cannot compute Inf-norm.

Implements Epetra_Operator.

Definition at line 209 of file Ifpack_ILU.h.

bool Ifpack_ILU::HasNormInf ( ) const
inlinevirtual

Returns false because this class cannot compute an Inf-norm.

Implements Epetra_Operator.

Definition at line 212 of file Ifpack_ILU.h.

bool Ifpack_ILU::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Epetra_Operator.

Definition at line 215 of file Ifpack_ILU.h.

const Epetra_Map& Ifpack_ILU::OperatorDomainMap ( ) const
inlinevirtual

Returns the Epetra_Map object associated with the domain of this operator.

Implements Epetra_Operator.

Definition at line 218 of file Ifpack_ILU.h.

const Epetra_Map& Ifpack_ILU::OperatorRangeMap ( ) const
inlinevirtual

Returns the Epetra_Map object associated with the range of this operator.

Implements Epetra_Operator.

Definition at line 221 of file Ifpack_ILU.h.

const Epetra_Comm& Ifpack_ILU::Comm ( ) const
inlinevirtual

Returns the Epetra_BlockMap object associated with the range of this matrix operator.

Implements Epetra_Operator.

Definition at line 224 of file Ifpack_ILU.h.

const Epetra_RowMatrix& Ifpack_ILU::Matrix ( ) const
inlinevirtual

Returns a reference to the matrix to be preconditioned.

Implements Ifpack_Preconditioner.

Definition at line 227 of file Ifpack_ILU.h.

std::ostream & Ifpack_ILU::Print ( std::ostream &  os) const
virtual

Prints on stream basic information about this object.

Implements Ifpack_Preconditioner.

Definition at line 634 of file Ifpack_ILU.cpp.

virtual int Ifpack_ILU::NumInitialize ( ) const
inlinevirtual

Returns the number of calls to Initialize().

Implements Ifpack_Preconditioner.

Definition at line 236 of file Ifpack_ILU.h.

virtual int Ifpack_ILU::NumCompute ( ) const
inlinevirtual

Returns the number of calls to Compute().

Implements Ifpack_Preconditioner.

Definition at line 242 of file Ifpack_ILU.h.

virtual int Ifpack_ILU::NumApplyInverse ( ) const
inlinevirtual

Returns the number of calls to ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 248 of file Ifpack_ILU.h.

virtual double Ifpack_ILU::InitializeTime ( ) const
inlinevirtual

Returns the time spent in Initialize().

Implements Ifpack_Preconditioner.

Definition at line 254 of file Ifpack_ILU.h.

virtual double Ifpack_ILU::ComputeTime ( ) const
inlinevirtual

Returns the time spent in Compute().

Implements Ifpack_Preconditioner.

Definition at line 260 of file Ifpack_ILU.h.

virtual double Ifpack_ILU::ApplyInverseTime ( ) const
inlinevirtual

Returns the time spent in ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 266 of file Ifpack_ILU.h.

virtual double Ifpack_ILU::InitializeFlops ( ) const
inlinevirtual

Returns the number of flops in the initialization phase.

Implements Ifpack_Preconditioner.

Definition at line 272 of file Ifpack_ILU.h.

virtual double Ifpack_ILU::ComputeFlops ( ) const
inlinevirtual

Returns the number of flops in the computation phase.

Implements Ifpack_Preconditioner.

Definition at line 277 of file Ifpack_ILU.h.

virtual double Ifpack_ILU::ApplyInverseFlops ( ) const
inlinevirtual

Returns the number of flops in the application of the preconditioner.

Implements Ifpack_Preconditioner.

Definition at line 282 of file Ifpack_ILU.h.

Ifpack_ILU& Ifpack_ILU::operator= ( const Ifpack_ILU )
inlineprivate

operator= (should never be used)

Definition at line 299 of file Ifpack_ILU.h.

void Ifpack_ILU::Destroy ( )
private

Destroys all internal data.

Definition at line 92 of file Ifpack_ILU.cpp.

int Ifpack_ILU::Solve ( bool  Trans,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
private

Returns the result of a Ifpack_ILU forward/back solve on a Epetra_MultiVector X in Y.

Parameters
InTrans -If true, solve transpose problem.
X- (In) A Epetra_MultiVector of dimension NumVectors to solve for.
OutY - (Out) A Epetra_MultiVector of dimension NumVectorscontaining result.
Returns
Integer error code, set to 0 if successful.

Definition at line 500 of file Ifpack_ILU.cpp.

int Ifpack_ILU::ComputeSetup ( )
private

Definition at line 115 of file Ifpack_ILU.cpp.

int Ifpack_ILU::InitAllValues ( const Epetra_RowMatrix A,
int  MaxNumEntries 
)
private
int Ifpack_ILU::LevelOfFill ( ) const
inlineprivate

Returns the level of fill.

Definition at line 324 of file Ifpack_ILU.h.

double Ifpack_ILU::RelaxValue ( ) const
inlineprivate

Get ILU(k) relaxation parameter.

Definition at line 327 of file Ifpack_ILU.h.

double Ifpack_ILU::AbsoluteThreshold ( ) const
inlineprivate

Get absolute threshold value.

Definition at line 330 of file Ifpack_ILU.h.

double Ifpack_ILU::RelativeThreshold ( ) const
inlineprivate

Get relative threshold value.

Definition at line 333 of file Ifpack_ILU.h.

int Ifpack_ILU::NumGlobalRows ( ) const
inlineprivate

Returns the number of global matrix rows.

Definition at line 337 of file Ifpack_ILU.h.

int Ifpack_ILU::NumGlobalCols ( ) const
inlineprivate

Returns the number of global matrix columns.

Definition at line 340 of file Ifpack_ILU.h.

int Ifpack_ILU::NumGlobalNonzeros ( ) const
inlineprivate

Returns the number of nonzero entries in the global graph.

Definition at line 343 of file Ifpack_ILU.h.

virtual int Ifpack_ILU::NumGlobalBlockDiagonals ( ) const
inlineprivatevirtual

Returns the number of diagonal entries found in the global input graph.

Definition at line 346 of file Ifpack_ILU.h.

long long Ifpack_ILU::NumGlobalRows64 ( ) const
inlineprivate

Definition at line 349 of file Ifpack_ILU.h.

long long Ifpack_ILU::NumGlobalCols64 ( ) const
inlineprivate

Definition at line 351 of file Ifpack_ILU.h.

long long Ifpack_ILU::NumGlobalNonzeros64 ( ) const
inlineprivate

Definition at line 353 of file Ifpack_ILU.h.

virtual long long Ifpack_ILU::NumGlobalBlockDiagonals64 ( ) const
inlineprivatevirtual

Definition at line 355 of file Ifpack_ILU.h.

int Ifpack_ILU::NumMyRows ( ) const
inlineprivate

Returns the number of local matrix rows.

Definition at line 358 of file Ifpack_ILU.h.

int Ifpack_ILU::NumMyCols ( ) const
inlineprivate

Returns the number of local matrix columns.

Definition at line 361 of file Ifpack_ILU.h.

int Ifpack_ILU::NumMyNonzeros ( ) const
inlineprivate

Returns the number of nonzero entries in the local graph.

Definition at line 364 of file Ifpack_ILU.h.

virtual int Ifpack_ILU::NumMyBlockDiagonals ( ) const
inlineprivatevirtual

Returns the number of diagonal entries found in the local input graph.

Definition at line 367 of file Ifpack_ILU.h.

virtual int Ifpack_ILU::NumMyDiagonals ( ) const
inlineprivatevirtual

Returns the number of nonzero diagonal values found in matrix.

Definition at line 370 of file Ifpack_ILU.h.

int Ifpack_ILU::IndexBase ( ) const
inlineprivate

Returns the index base for row and column indices for this graph.

Definition at line 374 of file Ifpack_ILU.h.

long long Ifpack_ILU::IndexBase64 ( ) const
inlineprivate

Definition at line 376 of file Ifpack_ILU.h.

const Ifpack_IlukGraph& Ifpack_ILU::Graph ( ) const
inlineprivate

Returns the address of the Ifpack_IlukGraph associated with this factored matrix.

Definition at line 379 of file Ifpack_ILU.h.

Epetra_RowMatrix& Ifpack_ILU::Matrix ( )
inlineprivate

Returns a reference to the matrix.

Definition at line 382 of file Ifpack_ILU.h.

Member Data Documentation

Teuchos::RefCountPtr<Epetra_RowMatrix> Ifpack_ILU::A_
private

Pointer to the Epetra_RowMatrix to factorize.

Definition at line 391 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Ifpack_IlukGraph> Ifpack_ILU::Graph_
private

Definition at line 392 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Epetra_CrsGraph> Ifpack_ILU::CrsGraph_
private

Definition at line 393 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Epetra_Map> Ifpack_ILU::IlukRowMap_
private

Definition at line 394 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Epetra_Map> Ifpack_ILU::IlukDomainMap_
private

Definition at line 395 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Epetra_Map> Ifpack_ILU::IlukRangeMap_
private

Definition at line 396 of file Ifpack_ILU.h.

const Epetra_Map* Ifpack_ILU::U_DomainMap_
private

Definition at line 397 of file Ifpack_ILU.h.

const Epetra_Map* Ifpack_ILU::L_RangeMap_
private

Definition at line 398 of file Ifpack_ILU.h.

const Epetra_Comm& Ifpack_ILU::Comm_
private

Definition at line 399 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Epetra_CrsMatrix> Ifpack_ILU::L_
private

Contains the L factors.

Definition at line 401 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Epetra_CrsMatrix> Ifpack_ILU::U_
private

Contains the U factors.

Definition at line 403 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Epetra_CrsGraph> Ifpack_ILU::L_Graph_
private

Definition at line 404 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Epetra_CrsGraph> Ifpack_ILU::U_Graph_
private

Definition at line 405 of file Ifpack_ILU.h.

Teuchos::RefCountPtr<Epetra_Vector> Ifpack_ILU::D_
private

Diagonal of factors.

Definition at line 407 of file Ifpack_ILU.h.

bool Ifpack_ILU::UseTranspose_
private

Definition at line 408 of file Ifpack_ILU.h.

int Ifpack_ILU::NumMyDiagonals_
private

Definition at line 410 of file Ifpack_ILU.h.

bool Ifpack_ILU::Allocated_
private

Definition at line 411 of file Ifpack_ILU.h.

bool Ifpack_ILU::ValuesInitialized_
private

Definition at line 412 of file Ifpack_ILU.h.

bool Ifpack_ILU::Factored_
private

Definition at line 413 of file Ifpack_ILU.h.

double Ifpack_ILU::RelaxValue_
private

Relaxation value.

Definition at line 415 of file Ifpack_ILU.h.

double Ifpack_ILU::Athresh_
private

absolute threshold

Definition at line 417 of file Ifpack_ILU.h.

double Ifpack_ILU::Rthresh_
private

relative threshold

Definition at line 419 of file Ifpack_ILU.h.

double Ifpack_ILU::Condest_
private

condition number estimate

Definition at line 421 of file Ifpack_ILU.h.

int Ifpack_ILU::LevelOfFill_
private

Level of fill.

Definition at line 423 of file Ifpack_ILU.h.

bool Ifpack_ILU::IsInitialized_
private

If true, the preconditioner has been successfully initialized.

Definition at line 425 of file Ifpack_ILU.h.

bool Ifpack_ILU::IsComputed_
private

If true, the preconditioner has been successfully computed.

Definition at line 427 of file Ifpack_ILU.h.

char Ifpack_ILU::Label_[160]
private

Label of this object.

Definition at line 429 of file Ifpack_ILU.h.

int Ifpack_ILU::NumInitialize_
private

Contains the number of successful calls to Initialize().

Definition at line 431 of file Ifpack_ILU.h.

int Ifpack_ILU::NumCompute_
private

Contains the number of successful call to Compute().

Definition at line 433 of file Ifpack_ILU.h.

int Ifpack_ILU::NumApplyInverse_
mutableprivate

Contains the number of successful call to ApplyInverse().

Definition at line 435 of file Ifpack_ILU.h.

double Ifpack_ILU::InitializeTime_
private

Contains the time for all successful calls to Initialize().

Definition at line 437 of file Ifpack_ILU.h.

double Ifpack_ILU::ComputeTime_
private

Contains the time for all successful calls to Compute().

Definition at line 439 of file Ifpack_ILU.h.

double Ifpack_ILU::ApplyInverseTime_
mutableprivate

Contains the time for all successful calls to ApplyInverse().

Definition at line 441 of file Ifpack_ILU.h.

double Ifpack_ILU::ComputeFlops_
private

Contains the number of flops for Compute().

Definition at line 443 of file Ifpack_ILU.h.

double Ifpack_ILU::ApplyInverseFlops_
mutableprivate

Contain sthe number of flops for ApplyInverse().

Definition at line 445 of file Ifpack_ILU.h.

Epetra_Time Ifpack_ILU::Time_
mutableprivate

Used for timing issues.

Definition at line 447 of file Ifpack_ILU.h.


The documentation for this class was generated from the following files: