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

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

#include <Ifpack_IC.h>

Inheritance diagram for Ifpack_IC:
Inheritance graph
[legend]

Public Member Functions

 Ifpack_IC (Epetra_RowMatrix *A)
 Ifpack_IC constuctor with variable number of indices per row. More...
 
virtual ~Ifpack_IC ()
 Ifpack_IC Destructor. More...
 
void SetAbsoluteThreshold (double Athresh)
 Set absolute threshold value. More...
 
void SetRelativeThreshold (double Rthresh)
 Set relative threshold value. More...
 
int SetParameters (Teuchos::ParameterList &parameterlis)
 Set parameters using a Teuchos::ParameterList object. More...
 
int SetParameter (const std::string, const int)
 
int SetParameter (const std::string, const double)
 
const Epetra_RowMatrixMatrix () const
 Returns a pointer to the matrix to be preconditioned. More...
 
Epetra_RowMatrixMatrix ()
 
bool IsInitialized () const
 Returns true if the preconditioner has been successfully initialized, false otherwise. 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...
 
int ComputeSetup ()
 
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_IC 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...
 
double GetAbsoluteThreshold ()
 Get absolute threshold value. More...
 
double GetRelativeThreshold ()
 Get relative threshold value. 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_VectorD () const
 Returns the address of the D factor associated with this factored matrix. More...
 
const Epetra_CrsMatrixU () const
 Returns the address of the U 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 the computation phase. More...
 
virtual double ApplyInverseFlops () const
 Returns the number of flops in the application of the preconditioner. More...
 

Private Member Functions

double LevelOfFill () const
 
double AbsoluteThreshold () const
 
double RelativeThreshold () const
 
double DropTolerance () const
 

Private Attributes

Teuchos::RefCountPtr
< Epetra_RowMatrix
A_
 
const Epetra_CommComm_
 
Teuchos::RefCountPtr
< Epetra_CrsMatrix
U_
 
Teuchos::RefCountPtr
< Epetra_Vector
D_
 
bool UseTranspose_
 
double Condest_
 
double Athresh_
 
double Rthresh_
 
double Droptol_
 
double Lfil_
 
void * Aict_
 
void * Lict_
 
double * Ldiag_
 
char Label_ [160]
 
bool IsInitialized_
 
bool IsComputed_
 
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...
 
Epetra_Time Time_
 Used for timing purposes. More...
 
double ComputeFlops_
 Contains the number of flops for Compute(). More...
 
double ApplyInverseFlops_
 Contain sthe number of flops for ApplyInverse(). 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...
 

Detailed Description

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

The Ifpack_IC class computes a threshold (not level) based incomplete

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

  1. Level of fill ratio. This defines the maximum number of entries per row/column in the factor, relative to the average number of nonzeros per row/col in A. The default value of 1.0 keeps the IC factor as sparse as A.

  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 80 of file Ifpack_IC.h.

Constructor & Destructor Documentation

Ifpack_IC::Ifpack_IC ( Epetra_RowMatrix A)

Ifpack_IC constuctor with variable number of indices per row.

Creates a Ifpack_IC object and allocates storage.

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

Definition at line 61 of file Ifpack_IC.cpp.

Ifpack_IC::~Ifpack_IC ( )
virtual

Ifpack_IC Destructor.

Definition at line 90 of file Ifpack_IC.cpp.

Member Function Documentation

void Ifpack_IC::SetAbsoluteThreshold ( double  Athresh)
inline

Set absolute threshold value.

Definition at line 97 of file Ifpack_IC.h.

void Ifpack_IC::SetRelativeThreshold ( double  Rthresh)
inline

Set relative threshold value.

Definition at line 100 of file Ifpack_IC.h.

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

Set parameters using a Teuchos::ParameterList object.

Implements Ifpack_Preconditioner.

Definition at line 110 of file Ifpack_IC.cpp.

int Ifpack_IC::SetParameter ( const std::string  ,
const int   
)
inline

Definition at line 112 of file Ifpack_IC.h.

int Ifpack_IC::SetParameter ( const std::string  ,
const double   
)
inline

Definition at line 116 of file Ifpack_IC.h.

const Epetra_RowMatrix& Ifpack_IC::Matrix ( ) const
inlinevirtual

Returns a pointer to the matrix to be preconditioned.

Implements Ifpack_Preconditioner.

Definition at line 121 of file Ifpack_IC.h.

Epetra_RowMatrix& Ifpack_IC::Matrix ( )
inline

Definition at line 126 of file Ifpack_IC.h.

bool Ifpack_IC::IsInitialized ( ) const
inlinevirtual

Returns true if the preconditioner has been successfully initialized, false otherwise.

Implements Ifpack_Preconditioner.

Definition at line 131 of file Ifpack_IC.h.

int Ifpack_IC::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 126 of file Ifpack_IC.cpp.

int Ifpack_IC::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 214 of file Ifpack_IC.cpp.

int Ifpack_IC::ComputeSetup ( )

Definition at line 139 of file Ifpack_IC.cpp.

bool Ifpack_IC::IsComputed ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 158 of file Ifpack_IC.h.

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

Returns the result of a Ifpack_IC 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 309 of file Ifpack_IC.cpp.

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

Implements Epetra_Operator.

Definition at line 350 of file Ifpack_IC.cpp.

double Ifpack_IC::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 370 of file Ifpack_IC.cpp.

double Ifpack_IC::Condest ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 191 of file Ifpack_IC.h.

double Ifpack_IC::GetAbsoluteThreshold ( )
inline

Get absolute threshold value.

Definition at line 199 of file Ifpack_IC.h.

double Ifpack_IC::GetRelativeThreshold ( )
inline

Get relative threshold value.

Definition at line 202 of file Ifpack_IC.h.

int Ifpack_IC::NumGlobalNonzeros ( ) const
inline

Returns the number of nonzero entries in the global graph.

Definition at line 206 of file Ifpack_IC.h.

long long Ifpack_IC::NumGlobalNonzeros64 ( ) const
inline

Definition at line 208 of file Ifpack_IC.h.

int Ifpack_IC::NumMyNonzeros ( ) const
inline

Returns the number of nonzero entries in the local graph.

Definition at line 211 of file Ifpack_IC.h.

const Epetra_Vector& Ifpack_IC::D ( ) const
inline

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

Definition at line 213 of file Ifpack_IC.h.

const Epetra_CrsMatrix& Ifpack_IC::U ( ) const
inline

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

Definition at line 216 of file Ifpack_IC.h.

int Ifpack_IC::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 230 of file Ifpack_IC.h.

double Ifpack_IC::NormInf ( ) const
inlinevirtual

Returns 0.0 because this class cannot compute Inf-norm.

Implements Epetra_Operator.

Definition at line 233 of file Ifpack_IC.h.

bool Ifpack_IC::HasNormInf ( ) const
inlinevirtual

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

Implements Epetra_Operator.

Definition at line 236 of file Ifpack_IC.h.

bool Ifpack_IC::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Epetra_Operator.

Definition at line 239 of file Ifpack_IC.h.

const Epetra_Map& Ifpack_IC::OperatorDomainMap ( ) const
inlinevirtual

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

Implements Epetra_Operator.

Definition at line 242 of file Ifpack_IC.h.

const Epetra_Map& Ifpack_IC::OperatorRangeMap ( ) const
inlinevirtual

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

Implements Epetra_Operator.

Definition at line 245 of file Ifpack_IC.h.

const Epetra_Comm& Ifpack_IC::Comm ( ) const
inlinevirtual

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

Implements Epetra_Operator.

Definition at line 248 of file Ifpack_IC.h.

const char* Ifpack_IC::Label ( ) const
inlinevirtual

Implements Epetra_Operator.

Definition at line 251 of file Ifpack_IC.h.

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

Definition at line 256 of file Ifpack_IC.h.

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

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

Implements Ifpack_Preconditioner.

Definition at line 385 of file Ifpack_IC.cpp.

virtual int Ifpack_IC::NumInitialize ( ) const
inlinevirtual

Returns the number of calls to Initialize().

Implements Ifpack_Preconditioner.

Definition at line 266 of file Ifpack_IC.h.

virtual int Ifpack_IC::NumCompute ( ) const
inlinevirtual

Returns the number of calls to Compute().

Implements Ifpack_Preconditioner.

Definition at line 272 of file Ifpack_IC.h.

virtual int Ifpack_IC::NumApplyInverse ( ) const
inlinevirtual

Returns the number of calls to ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 278 of file Ifpack_IC.h.

virtual double Ifpack_IC::InitializeTime ( ) const
inlinevirtual

Returns the time spent in Initialize().

Implements Ifpack_Preconditioner.

Definition at line 284 of file Ifpack_IC.h.

virtual double Ifpack_IC::ComputeTime ( ) const
inlinevirtual

Returns the time spent in Compute().

Implements Ifpack_Preconditioner.

Definition at line 290 of file Ifpack_IC.h.

virtual double Ifpack_IC::ApplyInverseTime ( ) const
inlinevirtual

Returns the time spent in ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 296 of file Ifpack_IC.h.

virtual double Ifpack_IC::InitializeFlops ( ) const
inlinevirtual

Returns the number of flops in the initialization phase.

Implements Ifpack_Preconditioner.

Definition at line 302 of file Ifpack_IC.h.

virtual double Ifpack_IC::ComputeFlops ( ) const
inlinevirtual

Returns the number of flops in the computation phase.

Implements Ifpack_Preconditioner.

Definition at line 307 of file Ifpack_IC.h.

virtual double Ifpack_IC::ApplyInverseFlops ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 312 of file Ifpack_IC.h.

double Ifpack_IC::LevelOfFill ( ) const
inlineprivate

Definition at line 320 of file Ifpack_IC.h.

double Ifpack_IC::AbsoluteThreshold ( ) const
inlineprivate

Definition at line 325 of file Ifpack_IC.h.

double Ifpack_IC::RelativeThreshold ( ) const
inlineprivate

Definition at line 330 of file Ifpack_IC.h.

double Ifpack_IC::DropTolerance ( ) const
inlineprivate

Definition at line 335 of file Ifpack_IC.h.

Member Data Documentation

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

Definition at line 340 of file Ifpack_IC.h.

const Epetra_Comm& Ifpack_IC::Comm_
private

Definition at line 341 of file Ifpack_IC.h.

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

Definition at line 342 of file Ifpack_IC.h.

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

Definition at line 343 of file Ifpack_IC.h.

bool Ifpack_IC::UseTranspose_
private

Definition at line 344 of file Ifpack_IC.h.

double Ifpack_IC::Condest_
private

Definition at line 346 of file Ifpack_IC.h.

double Ifpack_IC::Athresh_
private

Definition at line 347 of file Ifpack_IC.h.

double Ifpack_IC::Rthresh_
private

Definition at line 348 of file Ifpack_IC.h.

double Ifpack_IC::Droptol_
private

Definition at line 349 of file Ifpack_IC.h.

double Ifpack_IC::Lfil_
private

Definition at line 350 of file Ifpack_IC.h.

void* Ifpack_IC::Aict_
private

Definition at line 352 of file Ifpack_IC.h.

void* Ifpack_IC::Lict_
private

Definition at line 353 of file Ifpack_IC.h.

double* Ifpack_IC::Ldiag_
private

Definition at line 354 of file Ifpack_IC.h.

char Ifpack_IC::Label_[160]
private

Definition at line 355 of file Ifpack_IC.h.

bool Ifpack_IC::IsInitialized_
private

Definition at line 357 of file Ifpack_IC.h.

bool Ifpack_IC::IsComputed_
private

Definition at line 358 of file Ifpack_IC.h.

int Ifpack_IC::NumInitialize_
private

Contains the number of successful calls to Initialize().

Definition at line 361 of file Ifpack_IC.h.

int Ifpack_IC::NumCompute_
private

Contains the number of successful call to Compute().

Definition at line 363 of file Ifpack_IC.h.

int Ifpack_IC::NumApplyInverse_
mutableprivate

Contains the number of successful call to ApplyInverse().

Definition at line 365 of file Ifpack_IC.h.

double Ifpack_IC::InitializeTime_
private

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

Definition at line 368 of file Ifpack_IC.h.

double Ifpack_IC::ComputeTime_
private

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

Definition at line 370 of file Ifpack_IC.h.

double Ifpack_IC::ApplyInverseTime_
mutableprivate

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

Definition at line 372 of file Ifpack_IC.h.

Epetra_Time Ifpack_IC::Time_
mutableprivate

Used for timing purposes.

Definition at line 374 of file Ifpack_IC.h.

double Ifpack_IC::ComputeFlops_
private

Contains the number of flops for Compute().

Definition at line 377 of file Ifpack_IC.h.

double Ifpack_IC::ApplyInverseFlops_
mutableprivate

Contain sthe number of flops for ApplyInverse().

Definition at line 379 of file Ifpack_IC.h.


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