Ifpack Package Browser (Single Doxygen Collection)
Development
|
Ifpack_ILUT: A class for constructing and using an incomplete LU factorization of a given Epetra_RowMatrix. More...
#include <Ifpack_ILUT.h>
Ifpack_ILUT (const Epetra_RowMatrix *A) | |
Ifpack_ILUT constuctor with variable number of indices per row. More... | |
virtual | ~Ifpack_ILUT () |
Ifpack_ILUT Destructor. More... | |
int | SetParameters (Teuchos::ParameterList ¶meterlis) |
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_ILUT 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_Map & | OperatorDomainMap () const |
Returns the Epetra_Map object associated with the domain of this operator. More... | |
const Epetra_Map & | OperatorRangeMap () const |
Returns the Epetra_Map object associated with the range of this operator. More... | |
const Epetra_Comm & | Comm () const |
Returns the Epetra_BlockMap object associated with the range of this matrix operator. More... | |
const Epetra_RowMatrix & | Matrix () const |
Returns a reference to the matrix to be preconditioned. More... | |
const Epetra_CrsMatrix & | L () const |
Returns a reference to the L factor. More... | |
const Epetra_CrsMatrix & | U () 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_ILUT (const Ifpack_ILUT &RHS) | |
Copy constructor (should never be used) More... | |
Ifpack_ILUT & | operator= (const Ifpack_ILUT &) |
operator= (should never be used) More... | |
template<typename int_type > | |
int | TCompute () |
void | Destroy () |
Releases all allocated memory. More... | |
const Epetra_RowMatrix & | A_ |
reference to the matrix to be preconditioned. More... | |
const Epetra_Comm & | Comm_ |
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 | 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_Map > | SerialMap_ |
Additional Inherited Members |
Ifpack_ILUT: A class for constructing and using an incomplete LU factorization of a given Epetra_RowMatrix.
The Ifpack_ILUT class computes a "Relaxed" ILUT factorization with dual threshold dropping of small elements of a given Epetra_RowMatrix.
This implementation does not use the algorithm that is described in ifp_ilu. The algorithm drops entries in a row (i) of matrix A that are smaller than drop_tolerance even before the factorization of row i then computes the factorization for that row. This is different than the usual algorithm where the drop tolerance is applied to the factored rows.
The complete list of supported parameters is reported in page ifp_params.
Definition at line 81 of file Ifpack_ILUT.h.
Ifpack_ILUT::Ifpack_ILUT | ( | const Epetra_RowMatrix * | A | ) |
Ifpack_ILUT constuctor with variable number of indices per row.
Definition at line 66 of file Ifpack_ILUT.cpp.
|
virtual |
Ifpack_ILUT Destructor.
Definition at line 94 of file Ifpack_ILUT.cpp.
|
inlineprivate |
Copy constructor (should never be used)
Definition at line 319 of file Ifpack_ILUT.h.
|
virtual |
Set parameters using a Teuchos::ParameterList object.
Implements Ifpack_Preconditioner.
Definition at line 107 of file Ifpack_ILUT.cpp.
|
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.
In | A - User matrix to be factored. |
Implements Ifpack_Preconditioner.
Definition at line 144 of file Ifpack_ILUT.cpp.
|
inlinevirtual |
Returns true
if the preconditioner has been successfully initialized.
Implements Ifpack_Preconditioner.
Definition at line 113 of file Ifpack_ILUT.h.
|
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:
InitValues() must be called before the factorization can proceed.
Implements Ifpack_Preconditioner.
Definition at line 495 of file Ifpack_ILUT.cpp.
|
inlinevirtual |
If factor is completed, this query returns true, otherwise it returns false.
Implements Ifpack_Preconditioner.
Definition at line 130 of file Ifpack_ILUT.h.
|
virtual |
Returns the result of a Ifpack_ILUT forward/back solve on a Epetra_MultiVector X in Y.
X | - (In) A Epetra_MultiVector of dimension NumVectors to solve for. |
Y | - (Out) A Epetra_MultiVector of dimension NumVectorscontaining result. |
Implements Ifpack_Preconditioner.
Definition at line 537 of file Ifpack_ILUT.cpp.
|
virtual |
Implements Epetra_Operator.
Definition at line 584 of file Ifpack_ILUT.cpp.
|
virtual |
Computed the estimated condition number and returns the value.
Implements Ifpack_Preconditioner.
Definition at line 591 of file Ifpack_ILUT.cpp.
|
inlinevirtual |
Returns the computed estimated condition number, or -1.0 if no computed.
Implements Ifpack_Preconditioner.
Definition at line 154 of file Ifpack_ILUT.h.
|
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.
UseTranspose_in | - (In) If true, multiply by the transpose of operator, otherwise just use operator. |
Implements Epetra_Operator.
Definition at line 169 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns 0.0 because this class cannot compute Inf-norm.
Implements Epetra_Operator.
Definition at line 172 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns false because this class cannot compute an Inf-norm.
Implements Epetra_Operator.
Definition at line 175 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the current UseTranspose setting.
Implements Epetra_Operator.
Definition at line 178 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the Epetra_Map object associated with the domain of this operator.
Implements Epetra_Operator.
Definition at line 181 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the Epetra_Map object associated with the range of this operator.
Implements Epetra_Operator.
Definition at line 184 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
Implements Epetra_Operator.
Definition at line 187 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns a reference to the matrix to be preconditioned.
Implements Ifpack_Preconditioner.
Definition at line 190 of file Ifpack_ILUT.h.
|
inline |
Returns a reference to the L factor.
Definition at line 196 of file Ifpack_ILUT.h.
|
inline |
Returns a reference to the U factor.
Definition at line 199 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the label of this
object.
Implements Epetra_Operator.
Definition at line 202 of file Ifpack_ILUT.h.
|
inline |
Sets the label for this
object.
Definition at line 208 of file Ifpack_ILUT.h.
|
virtual |
Prints basic information on iostream. This function is used by operator<<.
Implements Ifpack_Preconditioner.
Definition at line 607 of file Ifpack_ILUT.cpp.
|
inlinevirtual |
Returns the number of calls to Initialize().
Implements Ifpack_Preconditioner.
Definition at line 218 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the number of calls to Compute().
Implements Ifpack_Preconditioner.
Definition at line 224 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the number of calls to ApplyInverse().
Implements Ifpack_Preconditioner.
Definition at line 230 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the time spent in Initialize().
Implements Ifpack_Preconditioner.
Definition at line 236 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the time spent in Compute().
Implements Ifpack_Preconditioner.
Definition at line 242 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the time spent in ApplyInverse().
Implements Ifpack_Preconditioner.
Definition at line 248 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the number of flops in the initialization phase.
Implements Ifpack_Preconditioner.
Definition at line 254 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the number of flops in the computation phase.
Implements Ifpack_Preconditioner.
Definition at line 259 of file Ifpack_ILUT.h.
|
inlinevirtual |
Returns the number of flops in the application of the preconditioner.
Implements Ifpack_Preconditioner.
Definition at line 264 of file Ifpack_ILUT.h.
|
inline |
Definition at line 269 of file Ifpack_ILUT.h.
|
inline |
Set relative threshold value.
Definition at line 274 of file Ifpack_ILUT.h.
|
inline |
Get absolute threshold value.
Definition at line 279 of file Ifpack_ILUT.h.
|
inline |
Get relative threshold value.
Definition at line 285 of file Ifpack_ILUT.h.
|
inline |
Gets the dropping tolerance.
Definition at line 291 of file Ifpack_ILUT.h.
|
inline |
Returns the number of nonzero entries in the global graph.
Definition at line 298 of file Ifpack_ILUT.h.
|
inline |
Definition at line 303 of file Ifpack_ILUT.h.
|
inline |
Returns the number of nonzero entries in the local graph.
Definition at line 309 of file Ifpack_ILUT.h.
|
inlineprivate |
operator= (should never be used)
Definition at line 326 of file Ifpack_ILUT.h.
|
private |
Definition at line 180 of file Ifpack_ILUT.cpp.
|
private |
Releases all allocated memory.
Definition at line 100 of file Ifpack_ILUT.cpp.
|
private |
reference to the matrix to be preconditioned.
Definition at line 341 of file Ifpack_ILUT.h.
|
private |
Reference to the communicator object.
Definition at line 343 of file Ifpack_ILUT.h.
|
private |
L factor.
Definition at line 345 of file Ifpack_ILUT.h.
|
private |
U factor.
Definition at line 347 of file Ifpack_ILUT.h.
|
private |
Condition number estimate.
Definition at line 349 of file Ifpack_ILUT.h.
|
private |
relaxation value
Definition at line 351 of file Ifpack_ILUT.h.
|
private |
Absolute threshold.
Definition at line 353 of file Ifpack_ILUT.h.
|
private |
Relative threshold.
Definition at line 355 of file Ifpack_ILUT.h.
|
private |
Level-of-fill.
Definition at line 357 of file Ifpack_ILUT.h.
|
private |
Discards all elements below this tolerance.
Definition at line 359 of file Ifpack_ILUT.h.
|
private |
Label for this
object.
Definition at line 361 of file Ifpack_ILUT.h.
|
private |
true
if this
object has been initialized
Definition at line 363 of file Ifpack_ILUT.h.
|
private |
true
if this
object has been computed
Definition at line 365 of file Ifpack_ILUT.h.
|
private |
true
if transpose has to be used.
Definition at line 367 of file Ifpack_ILUT.h.
|
private |
Number of local rows.
Definition at line 369 of file Ifpack_ILUT.h.
|
private |
Contains the number of successful calls to Initialize().
Definition at line 371 of file Ifpack_ILUT.h.
|
private |
Contains the number of successful call to Compute().
Definition at line 373 of file Ifpack_ILUT.h.
|
mutableprivate |
Contains the number of successful call to ApplyInverse().
Definition at line 375 of file Ifpack_ILUT.h.
|
private |
Contains the time for all successful calls to Initialize().
Definition at line 377 of file Ifpack_ILUT.h.
|
private |
Contains the time for all successful calls to Compute().
Definition at line 379 of file Ifpack_ILUT.h.
|
mutableprivate |
Contains the time for all successful calls to ApplyInverse().
Definition at line 381 of file Ifpack_ILUT.h.
|
private |
Contains the number of flops for Compute().
Definition at line 383 of file Ifpack_ILUT.h.
|
mutableprivate |
Contain sthe number of flops for ApplyInverse().
Definition at line 385 of file Ifpack_ILUT.h.
|
mutableprivate |
Used for timing purposed.
Definition at line 387 of file Ifpack_ILUT.h.
|
private |
Global number of nonzeros in L and U factors.
Definition at line 389 of file Ifpack_ILUT.h.
|
private |
Definition at line 390 of file Ifpack_ILUT.h.
|
private |
Definition at line 391 of file Ifpack_ILUT.h.