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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
virtual bool | destroyPreconditioner () const =0 |
Deletes the preconditioner. More... | |
virtual bool | recomputePreconditioner (const NOX::Epetra::Vector &x, Teuchos::ParameterList &linearSolverParams) const =0 |
Recalculates the preconditioner using an already allocated graph. More... | |
virtual PreconditionerReusePolicyType | getPreconditionerPolicy (bool advanceReuseCounter=true)=0 |
Evaluates the preconditioner policy at the current state. More... | |
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. More... | |
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. More... | |
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
where is the Jacobian, is the input vector, and is the result vector. Returns true if successful.
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
where is the Jacobian, is the input vector, and is the result vector.
The parameter list contains the linear solver options.
Implemented in NOX::Epetra::LinearSystemAztecOO.
Referenced by LOCA::BorderedSolver::EpetraAugmented::applyInverse().
|
pure virtual |
Applies Jacobian-Transpose to the given input vector and puts the answer in the result.
Computes
where is the Jacobian, is the input vector, and is the result vector. Returns true if successful.
Implemented in NOX::Epetra::LinearSystemAztecOO.
|
pure virtual |
Apply right preconditiong to the given input vector.
Let be a right preconditioner for the Jacobian ; in other words, is a matrix such that
Compute
where is the input vector and is the result vector.
If useTranspose is true, then the transpose of the preconditioner is applied:
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 |
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.
Referenced by LOCA::BorderedSolver::EpetraAugmented::applyInverse().
|
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.
Referenced by LOCA::Epetra::Group::applyComplexTransposeInverseMultiVector(), LOCA::BorderedSolver::EpetraAugmented::applyInverse(), LOCA::Epetra::Group::applyJacobianTransposeInverse(), and LOCA::Epetra::Group::applyJacobianTransposeInverseMultiVector().
|
pure virtual |
Return preconditioner operator.
Note: This should only be called if hasPreconditioner() returns true.
Implemented in NOX::Epetra::LinearSystemAztecOO.
Referenced by LOCA::BorderedSolver::EpetraAugmented::applyInverse().
|
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 |
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 preconditioner operator for solve.
Note: This should only be called if hasPreconditioner() returns true.
Implemented in NOX::Epetra::LinearSystemAztecOO.
Referenced by LOCA::BorderedSolver::EpetraAugmented::applyInverse().