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

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

#include <Ifpack_ICT.h>

Inheritance diagram for Ifpack_ICT:
Inheritance graph
[legend]

Public Member Functions

 Ifpack_ICT (const Epetra_RowMatrix *A)
 Ifpack_ICT constuctor with variable number of indices per row. More...
 
virtual ~Ifpack_ICT ()
 Ifpack_ICT Destructor. More...
 
int SetParameters (Teuchos::ParameterList &parameterlis)
 Set parameters using a Teuchos::ParameterList object. More...
 
const Epetra_RowMatrixMatrix () const
 Returns a reference to the matrix to be preconditioned. More...
 
bool IsInitialized () const
 Returns true is the preconditioner has been successfully initialized. More...
 
int Initialize ()
 Initialize L and U with values from user matrix A. 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_ICT 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)
 Returns the maximum over all the condition number estimate for each local ILU set of factors. More...
 
double Condest () const
 Returns the computed condition number estimate, or -1.0 if not computed. 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...
 
const Epetra_CrsMatrixH () const
 Returns the address of the D factor associated with this factored matrix. More...
 
const char * Label () const
 
int SetLabel (const char *Label_in)
 
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 all applications of Compute(). More...
 
virtual double ApplyInverseFlops () const
 Returns the number of flops in all applications of ApplyInverse(). More...
 
double LevelOfFill () const
 Returns the level-of-fill. More...
 
double AbsoluteThreshold () const
 Returns the absolute threshold. More...
 
double RelativeThreshold () const
 Returns the relative threshold. More...
 
double RelaxValue () const
 Returns the relaxation value. More...
 
double DropTolerance () const
 Returns the drop threshold. More...
 

Private Member Functions

 Ifpack_ICT (const Ifpack_ICT &rhs)
 Should not be used. More...
 
Ifpack_ICToperator= (const Ifpack_ICT &)
 Should not be used. More...
 
void Destroy ()
 Destroys all data associated to the preconditioner. More...
 
template<typename int_type >
int TCompute ()
 

Private Attributes

const Epetra_RowMatrixA_
 Reference to the matrix to be preconditioned, supposed symmetric. More...
 
const Epetra_CommComm_
 Reference to the communicator. More...
 
Teuchos::RefCountPtr
< Epetra_CrsMatrix
H_
 Contains the Cholesky factorization. More...
 
double Condest_
 Contains the estimate of the condition number, if -1.0 if not computed. More...
 
double Athresh_
 Absolute threshold. More...
 
double Rthresh_
 Relative threshold. More...
 
double LevelOfFill_
 Level of fill. More...
 
double DropTolerance_
 During factorization, drop all values below this. More...
 
double Relax_
 Relaxation value. More...
 
std::string Label_
 Label of this object. More...
 
bool IsInitialized_
 If true, the preconditioner has been successfully initialized. More...
 
bool IsComputed_
 If true, the preconditioner has been successfully computed. More...
 
bool UseTranspose_
 If true, use the transpose of the matrix. More...
 
int NumMyRows_
 Number of local rows in the matrix. 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 purposes. More...
 
long long GlobalNonzeros_
 Global number of nonzeros in L and U factors. More...
 
Teuchos::RefCountPtr
< Epetra_SerialComm
SerialComm_
 
Teuchos::RefCountPtr< Epetra_MapSerialMap_
 
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...
 

Detailed Description

Ifpack_ICT: A class for constructing and using an incomplete Cholesky factorization of a given Epetra_RowMatrix.

The Ifpack_ICT class computes a threshold based incomplete

LDL^T factorization of a given Epetra_RowMatrix. The factorization that is produced is a function of several parameters:

  1. Maximum number of entries per row/column in factor - The factorization will contain at most this number of nonzero terms in each row/column of the factorization.

  2. Diagonal perturbation - Prior to computing the factorization, it is possible to modify the diagonal entries of the matrix for which the factorization will be computing. If the absolute and relative perturbation values are zero and one, respectively, the factorization will be compute for the original user matrix A. Otherwise, the factorization will computed for a matrix that differs from the original user matrix in the diagonal values only. Details can be found in ifp_diag_pert.

Definition at line 83 of file Ifpack_ICT.h.

Constructor & Destructor Documentation

Ifpack_ICT::Ifpack_ICT ( const Epetra_RowMatrix A)

Ifpack_ICT constuctor with variable number of indices per row.

Creates a Ifpack_ICT object and allocates storage.

Parameters
InA - User matrix to be factored.
InGraph - Graph generated by Ifpack_IlukGraph.

Definition at line 63 of file Ifpack_ICT.cpp.

Ifpack_ICT::~Ifpack_ICT ( )
virtual

Ifpack_ICT Destructor.

Definition at line 91 of file Ifpack_ICT.cpp.

Ifpack_ICT::Ifpack_ICT ( const Ifpack_ICT rhs)
inlineprivate

Should not be used.

Definition at line 331 of file Ifpack_ICT.h.

Member Function Documentation

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

Set parameters using a Teuchos::ParameterList object.

Implements Ifpack_Preconditioner.

Definition at line 104 of file Ifpack_ICT.cpp.

const Epetra_RowMatrix& Ifpack_ICT::Matrix ( ) const
inlinevirtual

Returns a reference to the matrix to be preconditioned.

Implements Ifpack_Preconditioner.

Definition at line 110 of file Ifpack_ICT.h.

bool Ifpack_ICT::IsInitialized ( ) const
inlinevirtual

Returns true is the preconditioner has been successfully initialized.

Implements Ifpack_Preconditioner.

Definition at line 116 of file Ifpack_ICT.h.

int Ifpack_ICT::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 138 of file Ifpack_ICT.cpp.

int Ifpack_ICT::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 438 of file Ifpack_ICT.cpp.

bool Ifpack_ICT::IsComputed ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 142 of file Ifpack_ICT.h.

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

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

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

Implements Ifpack_Preconditioner.

Definition at line 455 of file Ifpack_ICT.cpp.

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

Implements Epetra_Operator.

Definition at line 494 of file Ifpack_ICT.cpp.

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

Returns the maximum over all the condition number estimate for each local ILU set of factors.

This functions computes a local condition number estimate on each processor and return the maximum over all processor of the estimate.

Parameters
InTrans -If true, solve transpose problem.
OutConditionNumberEstimate - The maximum across all processors of the infinity-norm estimate of the condition number of the inverse of LDU.

Implements Ifpack_Preconditioner.

Definition at line 502 of file Ifpack_ICT.cpp.

double Ifpack_ICT::Condest ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 175 of file Ifpack_ICT.h.

int Ifpack_ICT::NumGlobalNonzeros ( ) const
inline

Returns the number of nonzero entries in the global graph.

Definition at line 184 of file Ifpack_ICT.h.

long long Ifpack_ICT::NumGlobalNonzeros64 ( ) const
inline

Definition at line 186 of file Ifpack_ICT.h.

int Ifpack_ICT::NumMyNonzeros ( ) const
inline

Returns the number of nonzero entries in the local graph.

Definition at line 189 of file Ifpack_ICT.h.

const Epetra_CrsMatrix& Ifpack_ICT::H ( ) const
inline

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

Definition at line 192 of file Ifpack_ICT.h.

int Ifpack_ICT::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
InUseTranspose_in -If true, multiply by the transpose of operator, otherwise just use operator.
Returns
Always returns 0.

Implements Epetra_Operator.

Definition at line 206 of file Ifpack_ICT.h.

double Ifpack_ICT::NormInf ( ) const
inlinevirtual

Returns 0.0 because this class cannot compute Inf-norm.

Implements Epetra_Operator.

Definition at line 209 of file Ifpack_ICT.h.

bool Ifpack_ICT::HasNormInf ( ) const
inlinevirtual

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

Implements Epetra_Operator.

Definition at line 212 of file Ifpack_ICT.h.

bool Ifpack_ICT::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Epetra_Operator.

Definition at line 215 of file Ifpack_ICT.h.

const Epetra_Map& Ifpack_ICT::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_ICT.h.

const Epetra_Map& Ifpack_ICT::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_ICT.h.

const Epetra_Comm& Ifpack_ICT::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_ICT.h.

const char* Ifpack_ICT::Label ( ) const
inlinevirtual

Implements Epetra_Operator.

Definition at line 227 of file Ifpack_ICT.h.

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

Definition at line 232 of file Ifpack_ICT.h.

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

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

Implements Ifpack_Preconditioner.

Definition at line 518 of file Ifpack_ICT.cpp.

virtual int Ifpack_ICT::NumInitialize ( ) const
inlinevirtual

Returns the number of calls to Initialize().

Implements Ifpack_Preconditioner.

Definition at line 242 of file Ifpack_ICT.h.

virtual int Ifpack_ICT::NumCompute ( ) const
inlinevirtual

Returns the number of calls to Compute().

Implements Ifpack_Preconditioner.

Definition at line 248 of file Ifpack_ICT.h.

virtual int Ifpack_ICT::NumApplyInverse ( ) const
inlinevirtual

Returns the number of calls to ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 254 of file Ifpack_ICT.h.

virtual double Ifpack_ICT::InitializeTime ( ) const
inlinevirtual

Returns the time spent in Initialize().

Implements Ifpack_Preconditioner.

Definition at line 260 of file Ifpack_ICT.h.

virtual double Ifpack_ICT::ComputeTime ( ) const
inlinevirtual

Returns the time spent in Compute().

Implements Ifpack_Preconditioner.

Definition at line 266 of file Ifpack_ICT.h.

virtual double Ifpack_ICT::ApplyInverseTime ( ) const
inlinevirtual

Returns the time spent in ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 272 of file Ifpack_ICT.h.

virtual double Ifpack_ICT::InitializeFlops ( ) const
inlinevirtual

Returns the number of flops in the initialization phase.

Implements Ifpack_Preconditioner.

Definition at line 278 of file Ifpack_ICT.h.

virtual double Ifpack_ICT::ComputeFlops ( ) const
inlinevirtual

Returns the number of flops in all applications of Compute().

Implements Ifpack_Preconditioner.

Definition at line 284 of file Ifpack_ICT.h.

virtual double Ifpack_ICT::ApplyInverseFlops ( ) const
inlinevirtual

Returns the number of flops in all applications of ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 290 of file Ifpack_ICT.h.

double Ifpack_ICT::LevelOfFill ( ) const
inline

Returns the level-of-fill.

Note
: if 1.0, then the factored matrix contains approximatively the same number of elements of A.

Definition at line 299 of file Ifpack_ICT.h.

double Ifpack_ICT::AbsoluteThreshold ( ) const
inline

Returns the absolute threshold.

Definition at line 305 of file Ifpack_ICT.h.

double Ifpack_ICT::RelativeThreshold ( ) const
inline

Returns the relative threshold.

Definition at line 311 of file Ifpack_ICT.h.

double Ifpack_ICT::RelaxValue ( ) const
inline

Returns the relaxation value.

Definition at line 317 of file Ifpack_ICT.h.

double Ifpack_ICT::DropTolerance ( ) const
inline

Returns the drop threshold.

Definition at line 323 of file Ifpack_ICT.h.

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

Should not be used.

Definition at line 338 of file Ifpack_ICT.h.

void Ifpack_ICT::Destroy ( )
private

Destroys all data associated to the preconditioner.

Definition at line 97 of file Ifpack_ICT.cpp.

template<typename int_type >
int Ifpack_ICT::TCompute ( )
private

Definition at line 161 of file Ifpack_ICT.cpp.

Member Data Documentation

const Epetra_RowMatrix& Ifpack_ICT::A_
private

Reference to the matrix to be preconditioned, supposed symmetric.

Definition at line 347 of file Ifpack_ICT.h.

const Epetra_Comm& Ifpack_ICT::Comm_
private

Reference to the communicator.

Definition at line 349 of file Ifpack_ICT.h.

Teuchos::RefCountPtr<Epetra_CrsMatrix> Ifpack_ICT::H_
private

Contains the Cholesky factorization.

Definition at line 351 of file Ifpack_ICT.h.

double Ifpack_ICT::Condest_
private

Contains the estimate of the condition number, if -1.0 if not computed.

Definition at line 353 of file Ifpack_ICT.h.

double Ifpack_ICT::Athresh_
private

Absolute threshold.

Definition at line 355 of file Ifpack_ICT.h.

double Ifpack_ICT::Rthresh_
private

Relative threshold.

Definition at line 357 of file Ifpack_ICT.h.

double Ifpack_ICT::LevelOfFill_
private

Level of fill.

Definition at line 359 of file Ifpack_ICT.h.

double Ifpack_ICT::DropTolerance_
private

During factorization, drop all values below this.

Definition at line 361 of file Ifpack_ICT.h.

double Ifpack_ICT::Relax_
private

Relaxation value.

Definition at line 363 of file Ifpack_ICT.h.

std::string Ifpack_ICT::Label_
private

Label of this object.

Definition at line 365 of file Ifpack_ICT.h.

bool Ifpack_ICT::IsInitialized_
private

If true, the preconditioner has been successfully initialized.

Definition at line 367 of file Ifpack_ICT.h.

bool Ifpack_ICT::IsComputed_
private

If true, the preconditioner has been successfully computed.

Definition at line 369 of file Ifpack_ICT.h.

bool Ifpack_ICT::UseTranspose_
private

If true, use the transpose of the matrix.

Definition at line 371 of file Ifpack_ICT.h.

int Ifpack_ICT::NumMyRows_
private

Number of local rows in the matrix.

Definition at line 373 of file Ifpack_ICT.h.

int Ifpack_ICT::NumInitialize_
private

Contains the number of successful calls to Initialize().

Definition at line 375 of file Ifpack_ICT.h.

int Ifpack_ICT::NumCompute_
private

Contains the number of successful call to Compute().

Definition at line 377 of file Ifpack_ICT.h.

int Ifpack_ICT::NumApplyInverse_
mutableprivate

Contains the number of successful call to ApplyInverse().

Definition at line 379 of file Ifpack_ICT.h.

double Ifpack_ICT::InitializeTime_
private

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

Definition at line 381 of file Ifpack_ICT.h.

double Ifpack_ICT::ComputeTime_
private

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

Definition at line 383 of file Ifpack_ICT.h.

double Ifpack_ICT::ApplyInverseTime_
mutableprivate

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

Definition at line 385 of file Ifpack_ICT.h.

double Ifpack_ICT::ComputeFlops_
private

Contains the number of flops for Compute().

Definition at line 387 of file Ifpack_ICT.h.

double Ifpack_ICT::ApplyInverseFlops_
mutableprivate

Contain sthe number of flops for ApplyInverse().

Definition at line 389 of file Ifpack_ICT.h.

Epetra_Time Ifpack_ICT::Time_
mutableprivate

Used for timing purposes.

Definition at line 391 of file Ifpack_ICT.h.

long long Ifpack_ICT::GlobalNonzeros_
private

Global number of nonzeros in L and U factors.

Definition at line 393 of file Ifpack_ICT.h.

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

Definition at line 394 of file Ifpack_ICT.h.

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

Definition at line 395 of file Ifpack_ICT.h.


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