|
NOX Development
|
Pure virtual class interface for allowing different linear solvers to be used by the NOX::Epetra::Group. More...
#include <NOX_Epetra_LinearSystem.H>

Public Types | |
| enum | PreconditionerReusePolicyType { PRPT_REBUILD , PRPT_RECOMPUTE , PRPT_REUSE } |
| Determines handling of the preconditioner between nonlinear iterations. More... | |
Public Member Functions | |
| LinearSystem () | |
| Constructor. | |
| virtual | ~LinearSystem () |
| Destructor. | |
| virtual bool | applyJacobian (const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const =0 |
| Applies Jacobian to the given input vector and puts the answer in the result. | |
| virtual bool | applyJacobianTranspose (const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const =0 |
| Applies Jacobian-Transpose to the given input vector and puts the answer in the result. | |
| virtual bool | applyJacobianInverse (Teuchos::ParameterList ¶ms, const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result)=0 |
| Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result. | |
| virtual bool | applyRightPreconditioning (bool useTranspose, Teuchos::ParameterList ¶ms, const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const =0 |
| Apply right preconditiong to the given input vector. | |
| virtual Teuchos::RCP< NOX::Epetra::Scaling > | getScaling ()=0 |
| Get the scaling object. | |
| virtual void | resetScaling (const Teuchos::RCP< NOX::Epetra::Scaling > &s)=0 |
| Sets the diagonal scaling vector(s) used in scaling the linear system. | |
| virtual bool | computeJacobian (const NOX::Epetra::Vector &x)=0 |
| Evaluates the Jacobian based on the solution vector x. | |
| virtual bool | createPreconditioner (const NOX::Epetra::Vector &x, Teuchos::ParameterList &p, bool recomputeGraph) const =0 |
| Explicitly constructs a preconditioner based on the solution vector x and the parameter list p. | |
| virtual bool | destroyPreconditioner () const =0 |
| Deletes the preconditioner. | |
| virtual bool | recomputePreconditioner (const NOX::Epetra::Vector &x, Teuchos::ParameterList &linearSolverParams) const =0 |
| Recalculates the preconditioner using an already allocated graph. | |
| virtual PreconditionerReusePolicyType | getPreconditionerPolicy (bool advanceReuseCounter=true)=0 |
| Evaluates the preconditioner policy at the current state. | |
| virtual bool | isPreconditionerConstructed () const =0 |
| Indicates whether a preconditioner has been constructed. | |
| virtual bool | hasPreconditioner () const =0 |
| Indicates whether the linear system has a preconditioner. | |
| virtual Teuchos::RCP< const Epetra_Operator > | getJacobianOperator () const =0 |
| Return Jacobian operator. | |
| virtual Teuchos::RCP< Epetra_Operator > | getJacobianOperator ()=0 |
| Return Jacobian operator. | |
| virtual Teuchos::RCP< const Epetra_Operator > | getGeneratedPrecOperator () const =0 |
| Return preconditioner operator. | |
| virtual Teuchos::RCP< Epetra_Operator > | getGeneratedPrecOperator ()=0 |
| Return preconditioner operator. | |
| virtual void | setJacobianOperatorForSolve (const Teuchos::RCP< const Epetra_Operator > &solveJacOp)=0 |
| Set Jacobian operator for solve. | |
| virtual void | setPrecOperatorForSolve (const Teuchos::RCP< const Epetra_Operator > &solvePrecOp)=0 |
| Set preconditioner operator for solve. | |
| virtual int | getNumLinearSolves () |
| Statistics for number of times the linear solver has been called (def: 0) | |
| virtual int | getLinearItersLastSolve () |
| Statistics for number of iterations taken in last linear solve (def: 0) | |
| virtual int | getLinearItersTotal () |
| Statistics for cumulative number of iterations in all linear solve (def: 0) | |
| virtual double | getAchievedTol () |
| Statistics for the achieved tolerance of last linear solve (def: 0.0) | |
Pure virtual class interface for allowing different linear solvers to be used by the NOX::Epetra::Group.
Determines handling of the preconditioner between nonlinear iterations.
|
pure virtual |
Applies Jacobian to the given input vector and puts the answer in the result.
Computes
![\[ v = J u, \]](form_13.png)
where 


Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result.
Computes
![\[ v = J^{-1} u, \]](form_17.png)
where 


The parameter list contains the linear solver options.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Applies Jacobian-Transpose to the given input vector and puts the answer in the result.
Computes
![\[ v = J^T u, \]](form_16.png)
where 


Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Apply right preconditiong to the given input vector.
Let 


![\[ JM \approx I. \]](form_20.png)
Compute
![\[ u = M^{-1} v, \]](form_21.png)
where 

If useTranspose is true, then the transpose of the preconditioner is applied:
![\[ u = {M^{-1}}^T v, \]](form_22.png)
The transpose preconditioner is currently only required for Tensor methods.
The parameter list contains the linear solver options.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Evaluates the Jacobian based on the solution vector x.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Explicitly constructs a preconditioner based on the solution vector x and the parameter list p.
The user has the option of recomputing the graph when a new preconditioner is created. The NOX::Epetra::Group controls the isValid flag for the preconditioner and will control when to call this.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Deletes the preconditioner.
The NOX::Epetra::Group controls the isValid flag for the preconditioner and will control when to call this.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Return preconditioner operator.
Note: This should only be called if hasPreconditioner() returns true.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Return preconditioner operator.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Return Jacobian operator.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Return Jacobian operator.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Evaluates the preconditioner policy at the current state.
NOTE: This can change values between nonlienar iterations. It is not a static value.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Get the scaling object.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Indicates whether the linear system has a preconditioner.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Indicates whether a preconditioner has been constructed.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Recalculates the preconditioner using an already allocated graph.
Use this to compute a new preconditioner while using the same graph for the preconditioner. This avoids deleting and reallocating the memory required for the preconditioner and results in a big speed-up for large-scale jobs.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Sets the diagonal scaling vector(s) used in scaling the linear system.
See NOX::Epetra::Scaling for details on how to specify scaling of the linear system.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Set Jacobian operator for solve.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Set preconditioner operator for solve.
Note: This should only be called if hasPreconditioner() returns true.
Implemented in NOX::Epetra::LinearSystemAztecOO.