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_IKLU Class Reference

Ifpack_IKLU: A class for constructing and using an incomplete LU factorization of a given Epetra_RowMatrix. More...

#include <Ifpack_IKLU.h>

Inheritance diagram for Ifpack_IKLU:
Inheritance graph
[legend]
 Ifpack_IKLU (const Epetra_RowMatrix *A)
 Ifpack_IKLU constuctor with variable number of indices per row. More...
 
virtual ~Ifpack_IKLU ()
 Ifpack_IKLU Destructor. More...
 
int SetParameters (Teuchos::ParameterList &parameterlis)
 Set parameters using a Teuchos::ParameterList object. More...
 
int Initialize ()
 Initialize L and U with values from user matrix A. More...
 
bool IsInitialized () const
 Returns true if the preconditioner has been successfully initialized. More...
 
int Compute ()
 Compute IC factor 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 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. More...
 
int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 
double Condest (const Ifpack_CondestType CT=Ifpack_Cheap, const int MaxIters=1550, const double Tol=1e-9, Epetra_RowMatrix *Matrix_in=0)
 Computed the estimated condition number and returns the value. More...
 
double Condest () const
 Returns the computed estimated condition number, or -1.0 if no computed. More...
 
int SetUseTranspose (bool UseTranspose_in)
 If set true, transpose of this operator will be applied. 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...
 
const Epetra_CrsMatrixL () const
 Returns a reference to the L factor. More...
 
const Epetra_CrsMatrixU () const
 Returns a reference to the U factor. More...
 
const char * Label () const
 Returns the label of this object. More...
 
int SetLabel (const char *Label_in)
 Sets the label for this object. More...
 
virtual std::ostream & Print (std::ostream &os) const
 Prints basic information on iostream. This function is used by operator<<. 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...
 
double LevelOfFill () const
 
double RelaxValue () const
 Set relative threshold value. More...
 
double AbsoluteThreshold () const
 Get absolute threshold value. More...
 
double RelativeThreshold () const
 Get relative threshold value. More...
 
double DropTolerance () const
 Gets the dropping tolerance. More...
 
int NumGlobalNonzeros () const
 Returns the number of nonzero entries in the global graph. More...
 
long long NumGlobalNonzeros64 () const
 
int NumMyNonzeros () const
 Returns the number of nonzero entries in the local graph. More...
 
 Ifpack_IKLU (const Ifpack_IKLU &RHS)
 Copy constructor (should never be used) More...
 
Ifpack_IKLUoperator= (const Ifpack_IKLU &)
 operator= (should never be used) More...
 
void Destroy ()
 Releases all allocated memory. More...
 
const Epetra_RowMatrixA_
 reference to the matrix to be preconditioned. More...
 
const Epetra_CommComm_
 Reference to the communicator object. More...
 
Teuchos::RefCountPtr
< Epetra_CrsMatrix
L_
 L factor. More...
 
Teuchos::RefCountPtr
< Epetra_CrsMatrix
U_
 U factor. More...
 
double Condest_
 Condition number estimate. More...
 
double Relax_
 relaxation value More...
 
double Athresh_
 Absolute threshold. More...
 
double Rthresh_
 Relative threshold. More...
 
double LevelOfFill_
 Level-of-fill. More...
 
double DropTolerance_
 Discards all elements below this tolerance. More...
 
std::string Label_
 Label for this object. More...
 
bool IsInitialized_
 true if this object has been initialized More...
 
bool IsComputed_
 true if this object has been computed More...
 
bool UseTranspose_
 true if transpose has to be used. More...
 
int NumMyRows_
 Number of local rows. More...
 
int NumMyNonzeros_
 Number of local nonzeros. 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 purposed. More...
 
long long GlobalNonzeros_
 Global number of nonzeros in L and U factors. More...
 
Teuchos::RefCountPtr
< Epetra_SerialComm
SerialComm_
 
Teuchos::RefCountPtr< Epetra_MapSerialMap_
 
csrcsrA_
 Containers for the matrix storage and permutation. More...
 
csscssS_
 
csrncsrnN_
 Container for the L and U factor. More...
 

Additional Inherited Members

Detailed Description

Ifpack_IKLU: A class for constructing and using an incomplete LU factorization of a given Epetra_RowMatrix.

The Ifpack_IKLU class computes a "Relaxed" IKLU factorization with level k fill of a given Epetra_RowMatrix.

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
Heidi Thornquist, Org. 1437
Date
Last modified on 28-Nov-06.

Definition at line 79 of file Ifpack_IKLU.h.

Constructor & Destructor Documentation

Ifpack_IKLU::Ifpack_IKLU ( const Epetra_RowMatrix A)

Ifpack_IKLU constuctor with variable number of indices per row.

Definition at line 66 of file Ifpack_IKLU.cpp.

Ifpack_IKLU::~Ifpack_IKLU ( )
virtual

Ifpack_IKLU Destructor.

Definition at line 97 of file Ifpack_IKLU.cpp.

Ifpack_IKLU::Ifpack_IKLU ( const Ifpack_IKLU RHS)
inlineprivate

Copy constructor (should never be used)

Definition at line 317 of file Ifpack_IKLU.h.

Member Function Documentation

int Ifpack_IKLU::SetParameters ( Teuchos::ParameterList parameterlis)
virtual

Set parameters using a Teuchos::ParameterList object.

Implements Ifpack_Preconditioner.

Definition at line 116 of file Ifpack_IKLU.cpp.

int Ifpack_IKLU::Initialize ( )
virtual

Initialize L and U with values from user matrix A.

Copies values from the user's matrix into the nonzero pattern of L and U.

Parameters
InA - User matrix to be factored.
Warning
The graph of A must be identical to the graph passed in to Ifpack_IlukGraph constructor.

Implements Ifpack_Preconditioner.

Definition at line 153 of file Ifpack_IKLU.cpp.

bool Ifpack_IKLU::IsInitialized ( ) const
inlinevirtual

Returns true if the preconditioner has been successfully initialized.

Implements Ifpack_Preconditioner.

Definition at line 111 of file Ifpack_IKLU.h.

int Ifpack_IKLU::Compute ( )
virtual

Compute IC factor U using the specified graph, diagonal perturbation thresholds and relaxation parameters.

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

  1. Ifpack_IlukGraph specifying the structure of L and U.
  2. Value for the RILU(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 228 of file Ifpack_IKLU.cpp.

bool Ifpack_IKLU::IsComputed ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 128 of file Ifpack_IKLU.h.

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

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

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

Implements Ifpack_Preconditioner.

Definition at line 335 of file Ifpack_IKLU.cpp.

int Ifpack_IKLU::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Implements Epetra_Operator.

Definition at line 399 of file Ifpack_IKLU.cpp.

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

Computed the estimated condition number and returns the value.

Implements Ifpack_Preconditioner.

Definition at line 407 of file Ifpack_IKLU.cpp.

double Ifpack_IKLU::Condest ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 152 of file Ifpack_IKLU.h.

int Ifpack_IKLU::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 167 of file Ifpack_IKLU.h.

double Ifpack_IKLU::NormInf ( ) const
inlinevirtual

Returns 0.0 because this class cannot compute Inf-norm.

Implements Epetra_Operator.

Definition at line 170 of file Ifpack_IKLU.h.

bool Ifpack_IKLU::HasNormInf ( ) const
inlinevirtual

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

Implements Epetra_Operator.

Definition at line 173 of file Ifpack_IKLU.h.

bool Ifpack_IKLU::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Epetra_Operator.

Definition at line 176 of file Ifpack_IKLU.h.

const Epetra_Map& Ifpack_IKLU::OperatorDomainMap ( ) const
inlinevirtual

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

Implements Epetra_Operator.

Definition at line 179 of file Ifpack_IKLU.h.

const Epetra_Map& Ifpack_IKLU::OperatorRangeMap ( ) const
inlinevirtual

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

Implements Epetra_Operator.

Definition at line 182 of file Ifpack_IKLU.h.

const Epetra_Comm& Ifpack_IKLU::Comm ( ) const
inlinevirtual

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

Implements Epetra_Operator.

Definition at line 185 of file Ifpack_IKLU.h.

const Epetra_RowMatrix& Ifpack_IKLU::Matrix ( ) const
inlinevirtual

Returns a reference to the matrix to be preconditioned.

Implements Ifpack_Preconditioner.

Definition at line 188 of file Ifpack_IKLU.h.

const Epetra_CrsMatrix& Ifpack_IKLU::L ( ) const
inline

Returns a reference to the L factor.

Definition at line 194 of file Ifpack_IKLU.h.

const Epetra_CrsMatrix& Ifpack_IKLU::U ( ) const
inline

Returns a reference to the U factor.

Definition at line 197 of file Ifpack_IKLU.h.

const char* Ifpack_IKLU::Label ( ) const
inlinevirtual

Returns the label of this object.

Implements Epetra_Operator.

Definition at line 200 of file Ifpack_IKLU.h.

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

Sets the label for this object.

Definition at line 206 of file Ifpack_IKLU.h.

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

Prints basic information on iostream. This function is used by operator<<.

Implements Ifpack_Preconditioner.

Definition at line 423 of file Ifpack_IKLU.cpp.

virtual int Ifpack_IKLU::NumInitialize ( ) const
inlinevirtual

Returns the number of calls to Initialize().

Implements Ifpack_Preconditioner.

Definition at line 216 of file Ifpack_IKLU.h.

virtual int Ifpack_IKLU::NumCompute ( ) const
inlinevirtual

Returns the number of calls to Compute().

Implements Ifpack_Preconditioner.

Definition at line 222 of file Ifpack_IKLU.h.

virtual int Ifpack_IKLU::NumApplyInverse ( ) const
inlinevirtual

Returns the number of calls to ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 228 of file Ifpack_IKLU.h.

virtual double Ifpack_IKLU::InitializeTime ( ) const
inlinevirtual

Returns the time spent in Initialize().

Implements Ifpack_Preconditioner.

Definition at line 234 of file Ifpack_IKLU.h.

virtual double Ifpack_IKLU::ComputeTime ( ) const
inlinevirtual

Returns the time spent in Compute().

Implements Ifpack_Preconditioner.

Definition at line 240 of file Ifpack_IKLU.h.

virtual double Ifpack_IKLU::ApplyInverseTime ( ) const
inlinevirtual

Returns the time spent in ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 246 of file Ifpack_IKLU.h.

virtual double Ifpack_IKLU::InitializeFlops ( ) const
inlinevirtual

Returns the number of flops in the initialization phase.

Implements Ifpack_Preconditioner.

Definition at line 252 of file Ifpack_IKLU.h.

virtual double Ifpack_IKLU::ComputeFlops ( ) const
inlinevirtual

Returns the number of flops in the computation phase.

Implements Ifpack_Preconditioner.

Definition at line 257 of file Ifpack_IKLU.h.

virtual double Ifpack_IKLU::ApplyInverseFlops ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 262 of file Ifpack_IKLU.h.

double Ifpack_IKLU::LevelOfFill ( ) const
inline

Definition at line 267 of file Ifpack_IKLU.h.

double Ifpack_IKLU::RelaxValue ( ) const
inline

Set relative threshold value.

Definition at line 272 of file Ifpack_IKLU.h.

double Ifpack_IKLU::AbsoluteThreshold ( ) const
inline

Get absolute threshold value.

Definition at line 277 of file Ifpack_IKLU.h.

double Ifpack_IKLU::RelativeThreshold ( ) const
inline

Get relative threshold value.

Definition at line 283 of file Ifpack_IKLU.h.

double Ifpack_IKLU::DropTolerance ( ) const
inline

Gets the dropping tolerance.

Definition at line 289 of file Ifpack_IKLU.h.

int Ifpack_IKLU::NumGlobalNonzeros ( ) const
inline

Returns the number of nonzero entries in the global graph.

Definition at line 296 of file Ifpack_IKLU.h.

long long Ifpack_IKLU::NumGlobalNonzeros64 ( ) const
inline

Definition at line 301 of file Ifpack_IKLU.h.

int Ifpack_IKLU::NumMyNonzeros ( ) const
inline

Returns the number of nonzero entries in the local graph.

Definition at line 307 of file Ifpack_IKLU.h.

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

operator= (should never be used)

Definition at line 324 of file Ifpack_IKLU.h.

void Ifpack_IKLU::Destroy ( )
private

Releases all allocated memory.

Definition at line 103 of file Ifpack_IKLU.cpp.

Member Data Documentation

const Epetra_RowMatrix& Ifpack_IKLU::A_
private

reference to the matrix to be preconditioned.

Definition at line 336 of file Ifpack_IKLU.h.

const Epetra_Comm& Ifpack_IKLU::Comm_
private

Reference to the communicator object.

Definition at line 338 of file Ifpack_IKLU.h.

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

L factor.

Definition at line 340 of file Ifpack_IKLU.h.

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

U factor.

Definition at line 342 of file Ifpack_IKLU.h.

double Ifpack_IKLU::Condest_
private

Condition number estimate.

Definition at line 344 of file Ifpack_IKLU.h.

double Ifpack_IKLU::Relax_
private

relaxation value

Definition at line 346 of file Ifpack_IKLU.h.

double Ifpack_IKLU::Athresh_
private

Absolute threshold.

Definition at line 348 of file Ifpack_IKLU.h.

double Ifpack_IKLU::Rthresh_
private

Relative threshold.

Definition at line 350 of file Ifpack_IKLU.h.

double Ifpack_IKLU::LevelOfFill_
private

Level-of-fill.

Definition at line 352 of file Ifpack_IKLU.h.

double Ifpack_IKLU::DropTolerance_
private

Discards all elements below this tolerance.

Definition at line 354 of file Ifpack_IKLU.h.

std::string Ifpack_IKLU::Label_
private

Label for this object.

Definition at line 356 of file Ifpack_IKLU.h.

bool Ifpack_IKLU::IsInitialized_
private

true if this object has been initialized

Definition at line 358 of file Ifpack_IKLU.h.

bool Ifpack_IKLU::IsComputed_
private

true if this object has been computed

Definition at line 360 of file Ifpack_IKLU.h.

bool Ifpack_IKLU::UseTranspose_
private

true if transpose has to be used.

Definition at line 362 of file Ifpack_IKLU.h.

int Ifpack_IKLU::NumMyRows_
private

Number of local rows.

Definition at line 364 of file Ifpack_IKLU.h.

int Ifpack_IKLU::NumMyNonzeros_
private

Number of local nonzeros.

Definition at line 366 of file Ifpack_IKLU.h.

int Ifpack_IKLU::NumInitialize_
private

Contains the number of successful calls to Initialize().

Definition at line 368 of file Ifpack_IKLU.h.

int Ifpack_IKLU::NumCompute_
private

Contains the number of successful call to Compute().

Definition at line 370 of file Ifpack_IKLU.h.

int Ifpack_IKLU::NumApplyInverse_
mutableprivate

Contains the number of successful call to ApplyInverse().

Definition at line 372 of file Ifpack_IKLU.h.

double Ifpack_IKLU::InitializeTime_
private

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

Definition at line 374 of file Ifpack_IKLU.h.

double Ifpack_IKLU::ComputeTime_
private

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

Definition at line 376 of file Ifpack_IKLU.h.

double Ifpack_IKLU::ApplyInverseTime_
mutableprivate

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

Definition at line 378 of file Ifpack_IKLU.h.

double Ifpack_IKLU::ComputeFlops_
private

Contains the number of flops for Compute().

Definition at line 380 of file Ifpack_IKLU.h.

double Ifpack_IKLU::ApplyInverseFlops_
mutableprivate

Contain sthe number of flops for ApplyInverse().

Definition at line 382 of file Ifpack_IKLU.h.

Epetra_Time Ifpack_IKLU::Time_
mutableprivate

Used for timing purposed.

Definition at line 384 of file Ifpack_IKLU.h.

long long Ifpack_IKLU::GlobalNonzeros_
private

Global number of nonzeros in L and U factors.

Definition at line 386 of file Ifpack_IKLU.h.

Teuchos::RefCountPtr<Epetra_SerialComm> Ifpack_IKLU::SerialComm_
private

Definition at line 387 of file Ifpack_IKLU.h.

Teuchos::RefCountPtr<Epetra_Map> Ifpack_IKLU::SerialMap_
private

Definition at line 388 of file Ifpack_IKLU.h.

csr* Ifpack_IKLU::csrA_
private

Containers for the matrix storage and permutation.

Definition at line 391 of file Ifpack_IKLU.h.

css* Ifpack_IKLU::cssS_
private

Definition at line 392 of file Ifpack_IKLU.h.

csrn* Ifpack_IKLU::csrnN_
private

Container for the L and U factor.

Definition at line 394 of file Ifpack_IKLU.h.


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