46 #include "Ifpack_ConfigDefs.h" 47 #include "Ifpack_CondestType.h" 49 #include "Ifpack_Preconditioner.h" 50 #include "Epetra_Vector.h" 51 #include "Epetra_CrsMatrix.h" 52 #include "Epetra_Time.h" 53 #include "Teuchos_RefCountPtr.hpp" 55 class Epetra_RowMatrix;
56 class Epetra_SerialComm;
59 class Epetra_MultiVector;
107 int SetParameters(Teuchos::ParameterList& parameterlis);
118 return(IsInitialized_);
157 int ApplyInverse(
const Epetra_MultiVector& X, Epetra_MultiVector& Y)
const;
159 int Apply(
const Epetra_MultiVector& X, Epetra_MultiVector& Y)
const;
170 double Condest(
const Ifpack_CondestType CT = Ifpack_Cheap,
171 const int MaxIters = 1550,
172 const double Tol = 1e-9,
173 Epetra_RowMatrix* Matrix_in = 0);
183 #ifndef EPETRA_NO_32BIT_GLOBAL_INDICES 186 long long NumGlobalNonzeros64()
const {
return(H().NumGlobalNonzeros64());};
192 const Epetra_CrsMatrix&
H()
const {
return(*H_);};
206 int SetUseTranspose(
bool UseTranspose_in) {UseTranspose_ = UseTranspose_in;
return(0);};
224 const Epetra_Comm&
Comm()
const{
return(Comm_);};
227 const char* Label()
const 229 return(Label_.c_str());
232 int SetLabel(
const char* Label_in)
239 virtual std::ostream& Print(std::ostream& os)
const;
244 return(NumInitialize_);
256 return(NumApplyInverse_);
262 return(InitializeTime_);
268 return(ComputeTime_);
274 return(ApplyInverseTime_);
286 return(ComputeFlops_);
292 return(ApplyInverseFlops_);
301 return(LevelOfFill_);
325 return(DropTolerance_);
347 const Epetra_RowMatrix& A_;
349 const Epetra_Comm& Comm_;
351 Teuchos::RefCountPtr<Epetra_CrsMatrix> H_;
361 double DropTolerance_;
379 mutable int NumApplyInverse_;
381 double InitializeTime_;
385 mutable double ApplyInverseTime_;
387 double ComputeFlops_;
389 mutable double ApplyInverseFlops_;
391 mutable Epetra_Time Time_;
393 long long GlobalNonzeros_;
394 Teuchos::RefCountPtr<Epetra_SerialComm> SerialComm_;
395 Teuchos::RefCountPtr<Epetra_Map> SerialMap_;
397 template<
typename int_type>
double RelaxValue() const
Returns the relaxation value.
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
const Epetra_CrsMatrix & H() const
Returns the address of the D factor associated with this factored matrix.
double LevelOfFill() const
Returns the level-of-fill.
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
bool IsInitialized() const
Returns true is the preconditioner has been successfully initialized.
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
double DropTolerance() const
Returns the drop threshold.
Ifpack_ScalingType enumerable type.
virtual int NumInitialize() const
Returns the number of calls to Initialize().
virtual int NumCompute() const
Returns the number of calls to Compute().
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.
double Condest() const
Returns the computed condition number estimate, or -1.0 if not computed.
virtual double ComputeTime() const
Returns the time spent in Compute().
virtual double ApplyInverseTime() const
Returns the time spent in ApplyInverse().
virtual double ComputeFlops() const
Returns the number of flops in all applications of Compute().
bool IsComputed() const
If factor is completed, this query returns true, otherwise it returns false.
Ifpack_ICT: A class for constructing and using an incomplete Cholesky factorization of a given Epetra...
double RelativeThreshold() const
Returns the relative threshold.
bool UseTranspose() const
Returns the current UseTranspose setting.
const Epetra_RowMatrix & Matrix() const
Returns a reference to the matrix to be preconditioned.
int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied.
virtual double ApplyInverseFlops() const
Returns the number of flops in all applications of ApplyInverse().
virtual double InitializeTime() const
Returns the time spent in Initialize().
virtual double InitializeFlops() const
Returns the number of flops in the initialization phase.
virtual int NumApplyInverse() const
Returns the number of calls to ApplyInverse().
double AbsoluteThreshold() const
Returns the absolute threshold.