43 #ifndef IFPACK_POLYNOMIAL_H
44 #define IFPACK_POLYNOMIAL_H
46 #include "Ifpack_ConfigDefs.h"
47 #include "Ifpack_Preconditioner.h"
48 #include "Teuchos_RefCountPtr.hpp"
49 #include "Teuchos_LAPACK.hpp"
50 #include "Teuchos_SerialDenseMatrix.hpp"
65 #ifdef HAVE_IFPACK_EPETRAEXT
66 class EpetraExt_PointToBlockDiagPermute;
133 UseTranspose_ = UseTranspose_in;
174 virtual const char * Label()
const
176 return(Label_.c_str());
182 return(UseTranspose_);
204 return(IsInitialized_);
224 virtual double Condest(
const Ifpack_CondestType CT = Ifpack_Cheap,
225 const int MaxIters = 1550,
226 const double Tol = 1e-9,
239 virtual std::ostream&
Print(std::ostream & os)
const;
248 return(NumInitialize_);
260 return(NumApplyInverse_);
266 return(InitializeTime_);
272 return(ComputeTime_);
278 return(ApplyInverseTime_);
290 return(ComputeFlops_);
296 return(ApplyInverseFlops_);
305 const int MaximumIterations,
311 const int MaximumIterations,
312 double& lambda_min,
double& lambda_max);
314 #ifdef HAVE_IFPACK_EPETRAEXT
317 int CG(
const int MaximumIterations,
318 double& lambda_min,
double& lambda_max);
321 int PowerMethod(
const int MaximumIterations,
double& lambda_max);
327 const int MaximumIterations,
328 double& lambda_real_min,
double& lambda_real_max,
329 double& lambda_imag_min,
double& lambda_imag_max);
337 virtual void SetLabel();
363 mutable int NumApplyInverse_;
365 double InitializeTime_;
369 mutable double ApplyInverseTime_;
371 double ComputeFlops_;
373 mutable double ApplyInverseFlops_;
380 int LSPointsReal_, LSPointsImag_;
389 bool ComputeCondest_;
391 double RealEigRatio_, ImagEigRatio_;
399 double LambdaRealMin_, LambdaRealMax_, LambdaImagMin_, LambdaImagMax_;
401 double MinDiagonalValue_;
403 std::vector<double> coeff_;
411 long long NumGlobalRows_;
413 long long NumGlobalNonzeros_;
415 Teuchos::RefCountPtr<const Epetra_Operator> Operator_;
417 Teuchos::RefCountPtr<const Epetra_RowMatrix> Matrix_;
419 mutable Teuchos::RefCountPtr<Epetra_Vector> InvDiagonal_;
422 #ifdef HAVE_IFPACK_EPETRAEXT
423 Teuchos::ParameterList BlockList_;
425 Teuchos::RefCountPtr<EpetraExt_PointToBlockDiagPermute> InvBlockDiagonal_;
430 bool SolveNormalEquations_;
435 Teuchos::RefCountPtr<Epetra_Time> Time_;
437 bool ZeroStartingSolution_;
444 #endif // IFPACK_POLYNOMIAL_H
int GMRES(const Epetra_Operator &Operator, const Epetra_Vector &InvPointDiagonal, const int MaximumIterations, double &lambda_real_min, double &lambda_real_max, double &lambda_imag_min, double &lambda_imag_max)
Uses AztecOO's GMRES to estimate the height and width of the spectrum.
virtual bool HasNormInf() const
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
virtual int NumInitialize() const
Returns the number of calls to Initialize().
virtual double ComputeTime() const
Returns the time spent in Compute().
virtual bool UseTranspose() const
Returns the current UseTranspose setting.
virtual const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
virtual int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Applies the preconditioner to X, returns the result in Y.
virtual int Initialize()
Computes all it is necessary to initialize the preconditioner.
virtual int NumApplyInverse() const
Returns the number of calls to ApplyInverse().
virtual int SetUseTranspose(bool UseTranspose_in)
virtual double NormInf() const
Returns the infinity norm of the global matrix (not implemented)
static int CG(const Epetra_Operator &Operator, const Epetra_Vector &InvPointDiagonal, const int MaximumIterations, double &lambda_min, double &lambda_max)
Uses AztecOO's CG to estimate lambda_min and lambda_max.
virtual ~Ifpack_Polynomial()
Destructor.
virtual int Compute()
Computes the preconditioners.
virtual bool IsInitialized() const
Returns true if the preconditioner has been successfully initialized, false otherwise.
Ifpack_Preconditioner: basic class for preconditioning in Ifpack.
Ifpack_Polynomial(const Epetra_Operator *Matrix)
Ifpack_Polynomial constructor with given Epetra_Operator/Epetra_RowMatrix.
virtual int SetParameters(Teuchos::ParameterList &List)
Sets all the parameters for the preconditioner.
Ifpack_Polynomial: class for preconditioning with least squares polynomials in Ifpack.
virtual bool IsComputed() const
Returns true if the preconditioner has been successfully computed.
virtual const Epetra_RowMatrix & Matrix() const
Returns a pointer to the matrix to be preconditioned.
virtual std::ostream & Print(std::ostream &os) const
Prints object to an output stream.
virtual int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Applies the matrix to an Epetra_MultiVector.
virtual double InitializeFlops() const
Returns the number of flops in the initialization phase.
virtual int NumCompute() const
Returns the number of calls to Compute().
static int PowerMethod(const Epetra_Operator &Operator, const Epetra_Vector &InvPointDiagonal, const int MaximumIterations, double &LambdaMax)
Simple power method to compute lambda_max.
virtual double InitializeTime() const
Returns the time spent in Initialize().
virtual double ApplyInverseFlops() const
Returns the number of flops for the application of the preconditioner.
virtual double ApplyInverseTime() const
Returns the time spent in ApplyInverse().
virtual double ComputeFlops() const
Returns the number of flops in the computation phase.
virtual const Epetra_Comm & Comm() const
Returns a pointer to the Epetra_Comm communicator associated with this operator.
virtual const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
virtual double Condest() const
Returns the condition number estimate, or -1.0 if not computed.