|
Ifpack Package Browser (Single Doxygen Collection)
Development
|
Ifpack_Polynomial: class for preconditioning with least squares polynomials in Ifpack. More...
#include <Ifpack_Polynomial.h>

Public Member Functions | |
| virtual int | SetUseTranspose (bool UseTranspose_in) |
| Ifpack_Polynomial (const Epetra_Operator *Matrix) | |
| Ifpack_Polynomial constructor with given Epetra_Operator/Epetra_RowMatrix. More... | |
| Ifpack_Polynomial (const Epetra_RowMatrix *Matrix) | |
| Ifpack_Polynomial constructor with given Epetra_Operator/Epetra_RowMatrix. More... | |
| virtual | ~Ifpack_Polynomial () |
| Destructor. More... | |
| virtual int | Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
| Applies the matrix to an Epetra_MultiVector. More... | |
| virtual int | ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
| Applies the preconditioner to X, returns the result in Y. More... | |
| virtual double | NormInf () const |
| Returns the infinity norm of the global matrix (not implemented) More... | |
| virtual const char * | Label () const |
| virtual bool | UseTranspose () const |
| Returns the current UseTranspose setting. More... | |
| virtual bool | HasNormInf () const |
| Returns true if the this object can provide an approximate Inf-norm, false otherwise. More... | |
| virtual const Epetra_Comm & | Comm () const |
| Returns a pointer to the Epetra_Comm communicator associated with this operator. More... | |
| virtual const Epetra_Map & | OperatorDomainMap () const |
| Returns the Epetra_Map object associated with the domain of this operator. More... | |
| virtual const Epetra_Map & | OperatorRangeMap () const |
| Returns the Epetra_Map object associated with the range of this operator. More... | |
| virtual int | Initialize () |
| Computes all it is necessary to initialize the preconditioner. More... | |
| virtual bool | IsInitialized () const |
| Returns true if the preconditioner has been successfully initialized, false otherwise. More... | |
| virtual bool | IsComputed () const |
Returns true if the preconditioner has been successfully computed. More... | |
| virtual int | Compute () |
| Computes the preconditioners. More... | |
| virtual const Epetra_RowMatrix & | Matrix () const |
| Returns a pointer to the matrix to be preconditioned. More... | |
| virtual double | Condest (const Ifpack_CondestType CT=Ifpack_Cheap, const int MaxIters=1550, const double Tol=1e-9, Epetra_RowMatrix *Matrix_in=0) |
| Computes the condition number estimates and returns the value. More... | |
| virtual double | Condest () const |
| Returns the condition number estimate, or -1.0 if not computed. More... | |
| virtual int | SetParameters (Teuchos::ParameterList &List) |
| Sets all the parameters for the preconditioner. More... | |
| virtual std::ostream & | Print (std::ostream &os) const |
| Prints object to an output stream. 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 for the application of the preconditioner. More... | |
| 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. More... | |
| static int | PowerMethod (const Epetra_Operator &Operator, const Epetra_Vector &InvPointDiagonal, const int MaximumIterations, double &LambdaMax) |
| Simple power method to compute lambda_max. More... | |
| 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. More... | |
| bool | IsInitialized_ |
If true, the preconditioner has been computed successfully. More... | |
| bool | IsComputed_ |
If true, the preconditioner has been computed successfully. More... | |
| bool | IsIndefinite_ |
If true, have to compute polynomial for a spectrum with negative eigenvalues. More... | |
| bool | IsComplex_ |
If true, have to compute polynomial for a spectrum with nonzero imaginary part. 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... | |
| virtual void | SetLabel () |
| Sets the label. More... | |
| Ifpack_Polynomial (const Ifpack_Polynomial &) | |
| Copy constructor (PRIVATE, should not be used) More... | |
| Ifpack_Polynomial & | operator= (const Ifpack_Polynomial &) |
| operator = (PRIVATE, should not be used) More... | |
| int | PolyDegree_ |
| Contains the degree of the least squares polynomial. More... | |
| int | LSPointsReal_ |
| Contains the number of discretization points of the least squares problem. More... | |
| int | LSPointsImag_ |
| bool | UseTranspose_ |
If true, use the tranpose of Matrix_. More... | |
| double | Condest_ |
| Contains the estimated condition number. More... | |
| double | RealEigRatio_ |
| Contains the ratio such that the rectangular domain in the complex plane is [-LambdaRealMax_ / EigRatio_, LambdaRealMax_] x [-LambdaRealMax_ / ImagEigRatio_, LambdaRealMax_ / ImagEigRatio_]. More... | |
| double | ImagEigRatio_ |
| int | EigMaxIters_ |
| Max number of iterations to use in eigenvalue estimation (if automatic). More... | |
| std::string | Label_ |
| Contains the label of this object. More... | |
| double | LambdaRealMin_ |
| Bounds on the spectrum. More... | |
| double | LambdaRealMax_ |
| double | LambdaImagMin_ |
| double | LambdaImagMax_ |
| double | MinDiagonalValue_ |
| Contains the minimum value on the diagonal. More... | |
| std::vector< double > | coeff_ |
| coefficients of the polynomial More... | |
| int | NumMyRows_ |
| Number of local rows. More... | |
| int | NumMyNonzeros_ |
| Number of local nonzeros. More... | |
| long long | NumGlobalRows_ |
| Number of global rows. More... | |
| long long | NumGlobalNonzeros_ |
| Number of global nonzeros. More... | |
| Teuchos::RefCountPtr< const Epetra_Operator > | Operator_ |
| Pointers to the matrix to be preconditioned as an Epetra_Operator. More... | |
| Teuchos::RefCountPtr< const Epetra_RowMatrix > | Matrix_ |
| Pointers to the matrix to be preconditioned as an Epetra_RowMatrix. More... | |
| Teuchos::RefCountPtr < Epetra_Vector > | InvDiagonal_ |
Contains the inverse of diagonal elements of Matrix. More... | |
| bool | UseBlockMode_ |
| Use Block Preconditioning. More... | |
| bool | SolveNormalEquations_ |
| Run on the normal equations. More... | |
| bool | IsRowMatrix_ |
If true, the Operator_ is an Epetra_RowMatrix. More... | |
| Teuchos::RefCountPtr< Epetra_Time > | Time_ |
| Time object to track timing. More... | |
| bool | ZeroStartingSolution_ |
If true, the starting solution is always the zero vector. More... | |
Ifpack_Polynomial: class for preconditioning with least squares polynomials in Ifpack.
The Ifpack_Polynomial class enables the construction of preconditioners based on least squares polynomials for an Epetra_RowMatrix. Similar to Ifpack_Chebyshev, Ifpack_Polynomial is designed for indefinite linear systems.
The list of parameters is
Definition at line 99 of file Ifpack_Polynomial.h.
| Ifpack_Polynomial::Ifpack_Polynomial | ( | const Epetra_Operator * | Matrix | ) |
Ifpack_Polynomial constructor with given Epetra_Operator/Epetra_RowMatrix.
Creates an instance of Ifpack_Polynomial class.
| Matrix | - (In) Pointer to the operator to precondition. |
Definition at line 86 of file Ifpack_Polynomial.cpp.
| Ifpack_Polynomial::Ifpack_Polynomial | ( | const Epetra_RowMatrix * | Matrix | ) |
Ifpack_Polynomial constructor with given Epetra_Operator/Epetra_RowMatrix.
Creates an instance of Ifpack_Polynomial class.
| Matrix | - (In) Pointer to the matrix to precondition. |
Definition at line 135 of file Ifpack_Polynomial.cpp.
|
inlinevirtual |
Destructor.
Definition at line 121 of file Ifpack_Polynomial.h.
|
inlineprivate |
Copy constructor (PRIVATE, should not be used)
Definition at line 340 of file Ifpack_Polynomial.h.
|
inlinevirtual |
This flag can be used to apply the preconditioner to the transpose of the input operator.
Implements Epetra_Operator.
Definition at line 131 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Applies the matrix to an Epetra_MultiVector.
| X | - (In) A Epetra_MultiVector of dimension NumVectors to multiply with matrix. |
| Y | - (Out) A Epetra_MultiVector of dimension NumVectors containing the result. |
Implements Epetra_Operator.
Definition at line 252 of file Ifpack_Polynomial.cpp.
|
virtual |
Applies the preconditioner to X, returns the result in Y.
| X | - (In) A Epetra_MultiVector of dimension NumVectors to be preconditioned. |
| Y | - (InOut) A Epetra_MultiVector of dimension NumVectors containing result. |
Implements Ifpack_Preconditioner.
Definition at line 628 of file Ifpack_Polynomial.cpp.
|
inlinevirtual |
Returns the infinity norm of the global matrix (not implemented)
Implements Epetra_Operator.
Definition at line 166 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Implements Epetra_Operator.
Definition at line 174 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns the current UseTranspose setting.
Implements Epetra_Operator.
Definition at line 180 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns true if the this object can provide an approximate Inf-norm, false otherwise.
Implements Epetra_Operator.
Definition at line 186 of file Ifpack_Polynomial.h.
|
virtual |
Returns a pointer to the Epetra_Comm communicator associated with this operator.
Implements Epetra_Operator.
Definition at line 233 of file Ifpack_Polynomial.cpp.
|
virtual |
Returns the Epetra_Map object associated with the domain of this operator.
Implements Epetra_Operator.
Definition at line 239 of file Ifpack_Polynomial.cpp.
|
virtual |
Returns the Epetra_Map object associated with the range of this operator.
Implements Epetra_Operator.
Definition at line 245 of file Ifpack_Polynomial.cpp.
|
virtual |
Computes all it is necessary to initialize the preconditioner.
Implements Ifpack_Preconditioner.
Definition at line 273 of file Ifpack_Polynomial.cpp.
|
inlinevirtual |
Returns true if the preconditioner has been successfully initialized, false otherwise.
Implements Ifpack_Preconditioner.
Definition at line 202 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns true if the preconditioner has been successfully computed.
Implements Ifpack_Preconditioner.
Definition at line 208 of file Ifpack_Polynomial.h.
|
virtual |
Computes the preconditioners.
Implements Ifpack_Preconditioner.
Definition at line 307 of file Ifpack_Polynomial.cpp.
|
inlinevirtual |
Returns a pointer to the matrix to be preconditioned.
Implements Ifpack_Preconditioner.
Definition at line 218 of file Ifpack_Polynomial.h.
|
virtual |
Computes the condition number estimates and returns the value.
Implements Ifpack_Preconditioner.
Definition at line 606 of file Ifpack_Polynomial.cpp.
|
inlinevirtual |
Returns the condition number estimate, or -1.0 if not computed.
Implements Ifpack_Preconditioner.
Definition at line 230 of file Ifpack_Polynomial.h.
|
virtual |
Sets all the parameters for the preconditioner.
Implements Ifpack_Preconditioner.
Definition at line 177 of file Ifpack_Polynomial.cpp.
|
virtual |
Prints object to an output stream.
Implements Ifpack_Preconditioner.
Definition at line 548 of file Ifpack_Polynomial.cpp.
|
inlinevirtual |
Returns the number of calls to Initialize().
Implements Ifpack_Preconditioner.
Definition at line 246 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns the number of calls to Compute().
Implements Ifpack_Preconditioner.
Definition at line 252 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns the number of calls to ApplyInverse().
Implements Ifpack_Preconditioner.
Definition at line 258 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns the time spent in Initialize().
Implements Ifpack_Preconditioner.
Definition at line 264 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns the time spent in Compute().
Implements Ifpack_Preconditioner.
Definition at line 270 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns the time spent in ApplyInverse().
Implements Ifpack_Preconditioner.
Definition at line 276 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns the number of flops in the initialization phase.
Implements Ifpack_Preconditioner.
Definition at line 282 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns the number of flops in the computation phase.
Implements Ifpack_Preconditioner.
Definition at line 288 of file Ifpack_Polynomial.h.
|
inlinevirtual |
Returns the number of flops for the application of the preconditioner.
Implements Ifpack_Preconditioner.
Definition at line 294 of file Ifpack_Polynomial.h.
|
static |
Simple power method to compute lambda_max.
Definition at line 672 of file Ifpack_Polynomial.cpp.
|
static |
Uses AztecOO's CG to estimate lambda_min and lambda_max.
Definition at line 705 of file Ifpack_Polynomial.cpp.
| int Ifpack_Polynomial::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.
Definition at line 836 of file Ifpack_Polynomial.cpp.
|
privatevirtual |
Sets the label.
Definition at line 621 of file Ifpack_Polynomial.cpp.
|
inlineprivate |
operator = (PRIVATE, should not be used)
Definition at line 344 of file Ifpack_Polynomial.h.
|
private |
If true, the preconditioner has been computed successfully.
Definition at line 351 of file Ifpack_Polynomial.h.
|
private |
If true, the preconditioner has been computed successfully.
Definition at line 353 of file Ifpack_Polynomial.h.
|
private |
If true, have to compute polynomial for a spectrum with negative eigenvalues.
Definition at line 355 of file Ifpack_Polynomial.h.
|
private |
If true, have to compute polynomial for a spectrum with nonzero imaginary part.
Definition at line 357 of file Ifpack_Polynomial.h.
|
private |
Contains the number of successful calls to Initialize().
Definition at line 359 of file Ifpack_Polynomial.h.
|
private |
Contains the number of successful call to Compute().
Definition at line 361 of file Ifpack_Polynomial.h.
|
mutableprivate |
Contains the number of successful call to ApplyInverse().
Definition at line 363 of file Ifpack_Polynomial.h.
|
private |
Contains the time for all successful calls to Initialize().
Definition at line 365 of file Ifpack_Polynomial.h.
|
private |
Contains the time for all successful calls to Compute().
Definition at line 367 of file Ifpack_Polynomial.h.
|
mutableprivate |
Contains the time for all successful calls to ApplyInverse().
Definition at line 369 of file Ifpack_Polynomial.h.
|
private |
Contains the number of flops for Compute().
Definition at line 371 of file Ifpack_Polynomial.h.
|
mutableprivate |
Contain sthe number of flops for ApplyInverse().
Definition at line 373 of file Ifpack_Polynomial.h.
|
private |
Contains the degree of the least squares polynomial.
Definition at line 378 of file Ifpack_Polynomial.h.
|
private |
Contains the number of discretization points of the least squares problem.
Definition at line 380 of file Ifpack_Polynomial.h.
|
private |
Definition at line 380 of file Ifpack_Polynomial.h.
|
private |
If true, use the tranpose of Matrix_.
Definition at line 382 of file Ifpack_Polynomial.h.
|
private |
Contains the estimated condition number.
Definition at line 384 of file Ifpack_Polynomial.h.
|
private |
Contains the ratio such that the rectangular domain in the complex plane is [-LambdaRealMax_ / EigRatio_, LambdaRealMax_] x [-LambdaRealMax_ / ImagEigRatio_, LambdaRealMax_ / ImagEigRatio_].
Definition at line 393 of file Ifpack_Polynomial.h.
|
private |
Definition at line 393 of file Ifpack_Polynomial.h.
|
private |
Max number of iterations to use in eigenvalue estimation (if automatic).
Definition at line 395 of file Ifpack_Polynomial.h.
|
private |
Contains the label of this object.
Definition at line 397 of file Ifpack_Polynomial.h.
|
private |
Bounds on the spectrum.
Definition at line 399 of file Ifpack_Polynomial.h.
|
private |
Definition at line 399 of file Ifpack_Polynomial.h.
|
private |
Definition at line 399 of file Ifpack_Polynomial.h.
|
private |
Definition at line 399 of file Ifpack_Polynomial.h.
|
private |
Contains the minimum value on the diagonal.
Definition at line 401 of file Ifpack_Polynomial.h.
|
private |
coefficients of the polynomial
Definition at line 403 of file Ifpack_Polynomial.h.
|
private |
Number of local rows.
Definition at line 407 of file Ifpack_Polynomial.h.
|
private |
Number of local nonzeros.
Definition at line 409 of file Ifpack_Polynomial.h.
|
private |
Number of global rows.
Definition at line 411 of file Ifpack_Polynomial.h.
|
private |
Number of global nonzeros.
Definition at line 413 of file Ifpack_Polynomial.h.
|
private |
Pointers to the matrix to be preconditioned as an Epetra_Operator.
Definition at line 415 of file Ifpack_Polynomial.h.
|
private |
Pointers to the matrix to be preconditioned as an Epetra_RowMatrix.
Definition at line 417 of file Ifpack_Polynomial.h.
|
mutableprivate |
Contains the inverse of diagonal elements of Matrix.
Definition at line 419 of file Ifpack_Polynomial.h.
|
private |
Use Block Preconditioning.
Definition at line 421 of file Ifpack_Polynomial.h.
|
private |
Run on the normal equations.
Definition at line 430 of file Ifpack_Polynomial.h.
|
private |
If true, the Operator_ is an Epetra_RowMatrix.
Definition at line 433 of file Ifpack_Polynomial.h.
|
private |
Time object to track timing.
Definition at line 435 of file Ifpack_Polynomial.h.
|
private |
If true, the starting solution is always the zero vector.
Definition at line 437 of file Ifpack_Polynomial.h.
1.8.5