NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
NOX::Epetra::LinearSystemAztecOO Class Reference

Concrete implementation of NOX::Epetra::LinearSolver for AztecOO. More...

#include <NOX_Epetra_LinearSystem_AztecOO.H>

Inheritance diagram for NOX::Epetra::LinearSystemAztecOO:
Inheritance graph
[legend]
Collaboration diagram for NOX::Epetra::LinearSystemAztecOO:
Collaboration graph
[legend]

Public Member Functions

 LinearSystemAztecOO (Teuchos::ParameterList &printingParams, Teuchos::ParameterList &linearSolverParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &iReq, const NOX::Epetra::Vector &cloneVector, const Teuchos::RCP< NOX::Epetra::Scaling > scalingObject=Teuchos::null)
 Constructor with no Operators. More...
 
 LinearSystemAztecOO (Teuchos::ParameterList &printingParams, Teuchos::ParameterList &linearSolverParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &iReq, const Teuchos::RCP< NOX::Epetra::Interface::Jacobian > &iJac, const Teuchos::RCP< Epetra_Operator > &J, const NOX::Epetra::Vector &cloneVector, const Teuchos::RCP< NOX::Epetra::Scaling > scalingObject=Teuchos::null)
 Constructor with a user supplied Jacobian Operator only. More...
 
 LinearSystemAztecOO (Teuchos::ParameterList &printingParams, Teuchos::ParameterList &linearSolverParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &i, const Teuchos::RCP< NOX::Epetra::Interface::Preconditioner > &iPrec, const Teuchos::RCP< Epetra_Operator > &M, const NOX::Epetra::Vector &cloneVector, const Teuchos::RCP< NOX::Epetra::Scaling > scalingObject=Teuchos::null)
 Constructor with a user supplied Preconditioner Operator only. More...
 
 LinearSystemAztecOO (Teuchos::ParameterList &printingParams, Teuchos::ParameterList &linearSolverParams, const Teuchos::RCP< NOX::Epetra::Interface::Jacobian > &iJac, const Teuchos::RCP< Epetra_Operator > &J, const Teuchos::RCP< NOX::Epetra::Interface::Preconditioner > &iPrec, const Teuchos::RCP< Epetra_Operator > &M, const NOX::Epetra::Vector &cloneVector, const Teuchos::RCP< NOX::Epetra::Scaling > scalingObject=Teuchos::null)
 
virtual ~LinearSystemAztecOO ()
 Destructor.
 
virtual bool applyJacobian (const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const
 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
 Applies Jacobian-Transpose to the given input vector and puts the answer in the result. More...
 
virtual bool applyJacobianInverse (Teuchos::ParameterList &linearSolverParams, const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result)
 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 &linearSolverParams, const NOX::Epetra::Vector &input, NOX::Epetra::Vector &result) const
 Apply right preconditiong to the given input vector. More...
 
virtual bool createPreconditioner (const NOX::Epetra::Vector &x, Teuchos::ParameterList &linearSolverParams, bool recomputeGraph) const
 Explicitly constructs a preconditioner based on the solution vector x and the parameter list p. More...
 
virtual bool destroyPreconditioner () const
 Deletes all objects associated with the chosen preconditioner. This is called during linear solves and when the solution vector changes to reset the preconditioner.
 
virtual bool recomputePreconditioner (const NOX::Epetra::Vector &x, Teuchos::ParameterList &linearSolverParams) const
 Recalculates the preconditioner using an already allocated graph. More...
 
virtual
PreconditionerReusePolicyType 
getPreconditionerPolicy (bool advanceReuseCounter=true)
 Evaluates the preconditioner policy at the current state. More...
 
virtual void reset (Teuchos::ParameterList &linearSolverParams)
 Reset the linear solver parameters.
 
virtual Teuchos::RCP
< NOX::Epetra::Scaling
getScaling ()
 Get the scaling object.
 
void resetScaling (const Teuchos::RCP< NOX::Epetra::Scaling > &s)
 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.
 
virtual bool computeJacobian (const NOX::Epetra::Vector &x)
 Compute the Jacobian.
 
virtual Teuchos::RCP< const
NOX::Epetra::Interface::Jacobian
getJacobianInterface () const
 NOX::Interface::Jacobian accessor.
 
virtual Teuchos::RCP< const
NOX::Epetra::Interface::Preconditioner
getPrecInterface () const
 NOX::Interface::Preconditioiner accessor.
 
virtual bool isPreconditionerConstructed () const
 Indicates whether a preconditioner has been constructed.
 
virtual bool hasPreconditioner () const
 Indicates whether the linear system has a preconditioner.
 
virtual Teuchos::RCP< const
Epetra_Operator
getJacobianOperator () const
 Jacobian Epetra_Operator accessor.
 
virtual Teuchos::RCP
< Epetra_Operator
getJacobianOperator ()
 Jacobian Epetra_Operator accessor.
 
virtual Teuchos::RCP< const
Epetra_Operator
getPrecOperator () const
 Preconditioner Epetra_Operator accessor (only the base matrix if using an internal preconditioner - aztecoo or ifpack).
 
virtual Teuchos::RCP< const
Epetra_Operator
getGeneratedPrecOperator () const
 Return preconditioner operator generated and stored in AztecOO. More...
 
virtual Teuchos::RCP
< Epetra_Operator
getGeneratedPrecOperator ()
 Return preconditioner operator generated and stored in AztecOO.
 
double getTimeCreatePreconditioner () const
 Returns the total time (sec.) spent in createPreconditioner().
 
double getTimeApplyJacobianInverse () const
 Returns the total time (sec.) spent in applyJacobianInverse().
 
virtual void setJacobianOperatorForSolve (const Teuchos::RCP< const Epetra_Operator > &solveJacOp)
 Set Jacobian operator for solve.
 
virtual void setPrecOperatorForSolve (const Teuchos::RCP< const Epetra_Operator > &solvePrecOp)
 Set preconditioner operator for solve. More...
 
- Public Member Functions inherited from NOX::Epetra::LinearSystem
 LinearSystem ()
 Constructor.
 
virtual ~LinearSystem ()
 Destructor.
 
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)
 

Protected Types

enum  OperatorType { EpetraOperator, EpetraRowMatrix, EpetraVbrMatrix, EpetraCrsMatrix }
 List of types of epetra objects that can be used for the Jacobian and/or Preconditioner. More...
 
enum  PreconditionerMatrixSourceType { UseJacobian, SeparateMatrix }
 Source of the RowMatrix if using an AztecOO native preconditioner.
 
enum  PreconditionerType {
  None_, AztecOO_, Ifpack_, NewIfpack_,
  ML_, UserDefined_
}
 

Protected Member Functions

virtual void setAztecOptions (Teuchos::ParameterList &lsParams, AztecOO &aztec) const
 Parse the parameter list and set the corresponding options in the AztecOO solver onject.
 
virtual bool createJacobianOperator (Teuchos::ParameterList &printParams, Teuchos::ParameterList &lsParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &iReq, const NOX::Epetra::Vector &cloneVector)
 Creates an internally owned Epetra_Operator for the Jacobian.
 
virtual bool createPrecOperator (Teuchos::ParameterList &printParams, Teuchos::ParameterList &lsParams, const Teuchos::RCP< NOX::Epetra::Interface::Required > &iReq, const NOX::Epetra::Vector &cloneVector)
 Creates an internally owned Epetra_Operator for the Preconditioner.
 
virtual bool checkPreconditionerValidity ()
 Checks to make sure that the supplied operators are valid for the requested preconditioning options set in the parameter list.
 
virtual bool createIfpackPreconditioner (Teuchos::ParameterList &p) const
 Allocates the objects required for using ifpack preconditioners (NOX::Epetra::Group::ifpackGraph and NOX::Epetra::Group::ifpackPreconditioner). This is called from NOX::Epetra::Group::computePreconditioner().
 
virtual bool createNewIfpackPreconditioner (Teuchos::ParameterList &p) const
 Allocates the objects required for using the new version of ifpack preconditioners via the Ifpack Create factory. This is called from NOX::Epetra::Group::computePreconditioner().
 
virtual OperatorType getOperatorType (const Epetra_Operator &o)
 Deletes the AztecOO solver object. This is called when the solution vector for the group is changed. The preconditioner is no longer valid so the solver and preconditioner are destroyed by a call to this method. More...
 
virtual void setAztecOOJacobian () const
 Sets the epetra Jacobian operator in the AztecOO object. More...
 
virtual void setAztecOOPreconditioner () const
 Sets the epetra Preconditioner operator in the AztecOO object. More...
 
virtual void throwError (const std::string &functionName, const std::string &errorMsg) const
 
virtual void precError (int error_code, const std::string &nox_function, const std::string &prec_type, const std::string &prec_function) const
 Prints a warning for ifpack preconditioner failures (error_code != 0).
 

Protected Attributes

NOX::Utils utils
 Printing Utilities object.
 
Teuchos::RCP
< NOX::Epetra::Interface::Jacobian
jacInterfacePtr
 Reference to the user supplied Jacobian interface functions.
 
OperatorType jacType
 Type of operator for the Jacobian.
 
Teuchos::RCP< Epetra_OperatorjacPtr
 Pointer to the Jacobian operator.
 
Teuchos::RCP
< NOX::Epetra::Interface::Preconditioner
precInterfacePtr
 Reference to the user supplied preconditioner interface functions.
 
OperatorType precType
 Type of operator for the preconditioner.
 
Teuchos::RCP< Epetra_OperatorprecPtr
 Pointer to the preconditioner operator.
 
PreconditionerMatrixSourceType precMatrixSource
 
PreconditionerType precAlgorithm
 
Teuchos::RCP< AztecOO > aztecSolverPtr
 Aztec solver object.
 
Teuchos::RCP< Ifpack_IlukGraph > ifpackGraphPtr
 Stores the ifpack graph. Mutable since the applyRightPreconditioner() is a const method.
 
Teuchos::RCP< Ifpack_CrsRiluk > ifpackPreconditionerPtr
 Stores the ifpack preconditioner. Mutable since the applyRightPreconditioner() is a const method.
 
Teuchos::RCP
< Ifpack_Preconditioner > 
newIfpackPreconditionerPtr
 Stores the new ifpack preconditioner. Mutable since the applyRightPreconditioner() is a const method.
 
Teuchos::RCP
< NOX::Epetra::Scaling
scaling
 Scaling object supplied by the user.
 
Teuchos::RCP< NOX::Epetra::VectortmpVectorPtr
 An extra temporary vector, only allocated if needed.
 
double conditionNumberEstimate
 
bool isPrecConstructed
 True if the preconditioner has been computed.
 
bool outputSolveDetails
 If set to true, solver information is printed to the "Output" sublist of the "Linear Solver" list.
 
bool zeroInitialGuess
 Zero out the initial guess for linear solves performed through applyJacobianInverse calls (i.e. zero out the result vector before the linear solve).
 
bool manualScaling
 Stores the parameter "Compute Scaling Manually".
 
PreconditionerReusePolicyType precReusePolicy
 Policy for how to handle the preconditioner between nonlineaer iterations.
 
int precQueryCounter
 Counter for number of times called since reset or construction.
 
int maxAgeOfPrec
 Parameter to determine whether or not to recompute Preconditioner.
 
Epetra_Time timer
 Epetra_Time object.
 
double timeCreatePreconditioner
 Total time spent in createPreconditioner (sec.).
 
double timeApplyJacbianInverse
 Total time spent in applyJacobianInverse (sec.).
 
Teuchos::RCP< Epetra_OperatorsolvePrecOpPtr
 Preconditioner operator that will be used in solves.
 
bool throwErrorOnPrecFailure
 If true, any preconditioner error will cause a throw instead of a warning.
 

Additional Inherited Members

- Public Types inherited from NOX::Epetra::LinearSystem
enum  PreconditionerReusePolicyType { PRPT_REBUILD, PRPT_RECOMPUTE, PRPT_REUSE }
 Determines handling of the preconditioner between nonlinear iterations. More...
 

Detailed Description

Concrete implementation of NOX::Epetra::LinearSolver for AztecOO.

This solver provides the linear algebra services provided through the AztecOO parallel iterative linear solver.

The NOX::Epetra::LinearSystemAztecOO object provides a flexible and efficient way to interface an Epetra based application code to the Aztec linear solver. This class handles construction of both the preconditioners and AztecOO solver. All options are determined through parameter lists and the basic constructors.

Constructing a Linear System

There are four different constructors that can be used. The difference between constructors is based on whether the user supplies a Jacobian, a preconditioner, neither or both.

If a Jacobian is not supplied then this object can create an internally constructed Jacobian based on a Finite Difference or Matrif-Free object. The user can specify which type of object to use by setting the parameter "Jacobian Operator" in the parameter list. The choices are "Matrix-Free" or "Finite Difference".

The user can supply their own preconditioner as an Epetra_Operator, or they can supply their own matrix (an Epetra_RowMatrix derived object) that can be used by one of the internal preconditioner libraries (currently aztecoo or ifpack). If they supply their own preconditioner the object must implement the Epetra_Operator::ApplyInverse method. This is the method called during the linear solve to introduce preconditoning into aztecoo. If the user supplies a matrix to be used with an internal preconditioner, it must be derived from the Epetra_RowMatrix class and must implement all functionality in the Epetra_RowMatrix. If a Preconditioner is not supplied, then this object can create an internal preconditioner matrix by finite differencing or it can use the Jacobian operator if the Jacobian derives from the Epetra_RowMatrix class. The user can specify which type of object to use by setting the parameter "Preconditioner Operator" in the parameter list. The choices are "Use Jacobian" or "Finite Difference".

The Jacobian and preconditioner each require an interface to update the state of the operator with respect to the solution vector and any other parameters. There are three interfaces that can be implemented, NOX::Epetra::Interface::Required, NOX::Epetra::Interface::Jacobian, and NOX::Epetra::Interface::Preconditioner.

NOX::Epetra::Interface::Required supplies the computeF() function so codes can tell NOX what the nonlinear equations are. This is the minimum requirement to run nox through the epetra interface. LinearSolverAztecOO requires this in some constructors so that if a Jacobian or preconditoner is not supplied, it will use computeF from the Required interface to estimate the Jacobian or preconditioner via finite differences or directional derivatives.

NOX::Epetra::Interface::Jacobian is used for updating a user supplied Jacobian opertor with respect to the solution vector and any other parameters. It is required only in constructors in which a user supplies a Jacobian operator.

NOX::Epetra::Interface::Preconditioner is used for updating a user supplied preconditioner opertor/matrix with respect to the solution vector and any other parameters. It is required only in constructors in which a user supplies a preconditioner operator.

"Linear Solver" sublist parameters

A Teuchos::ParameterList called linearSolverParams is required in the various constructors and during some method calls such as applyJacobianInverse() and applyRightPreconditioning(). Typically, this list is the "Linear Solver" sublist found in the nox parameter list. The following parameters can be set in the linear solver sublist and are valid for the NOX::Epetra::LinearSolverAztecOO object:

"Output" sublist

The parameter list passed in during calls to ApplyJacobianInverse() will have an "Output" sublist created that contains the following parameters if the flag "Output Solver Details" is set to true:

Member Enumeration Documentation

List of types of epetra objects that can be used for the Jacobian and/or Preconditioner.

Enumerator
EpetraOperator 

An Epetra_Operator derived object.

EpetraRowMatrix 

An Epetra_RowMatrix derived object.

EpetraVbrMatrix 

An Epetra_VbrMatrix object.

EpetraCrsMatrix 

An Epetra_CrsMatrix object.

Constructor & Destructor Documentation

NOX::Epetra::LinearSystemAztecOO::LinearSystemAztecOO ( Teuchos::ParameterList printingParams,
Teuchos::ParameterList linearSolverParams,
const Teuchos::RCP< NOX::Epetra::Interface::Required > &  iReq,
const NOX::Epetra::Vector cloneVector,
const Teuchos::RCP< NOX::Epetra::Scaling scalingObject = Teuchos::null 
)

Constructor with no Operators.

Jacobian Operator will be constructed internally based on the parameter "Jacobian Operator". Defaults to using a NOX::Epetra::MatrixFree object.

References aztecSolverPtr, createJacobianOperator(), createPrecOperator(), Teuchos::rcp(), reset(), and tmpVectorPtr.

NOX::Epetra::LinearSystemAztecOO::LinearSystemAztecOO ( Teuchos::ParameterList printingParams,
Teuchos::ParameterList linearSolverParams,
const Teuchos::RCP< NOX::Epetra::Interface::Required > &  iReq,
const Teuchos::RCP< NOX::Epetra::Interface::Jacobian > &  iJac,
const Teuchos::RCP< Epetra_Operator > &  J,
const NOX::Epetra::Vector cloneVector,
const Teuchos::RCP< NOX::Epetra::Scaling scalingObject = Teuchos::null 
)

Constructor with a user supplied Jacobian Operator only.

Either there is no preconditioning or the preconditioner will be used/created internally. The Jacobian (if derived from an Epetra_RowMatrix class can be used with an internal preconditioner. See the parameter key "Preconditioner Operator" for more details.

References aztecSolverPtr, createPrecOperator(), getOperatorType(), jacPtr, jacType, Teuchos::rcp(), reset(), and tmpVectorPtr.

NOX::Epetra::LinearSystemAztecOO::LinearSystemAztecOO ( Teuchos::ParameterList printingParams,
Teuchos::ParameterList linearSolverParams,
const Teuchos::RCP< NOX::Epetra::Interface::Required > &  i,
const Teuchos::RCP< NOX::Epetra::Interface::Preconditioner > &  iPrec,
const Teuchos::RCP< Epetra_Operator > &  M,
const NOX::Epetra::Vector cloneVector,
const Teuchos::RCP< NOX::Epetra::Scaling scalingObject = Teuchos::null 
)

Constructor with a user supplied Preconditioner Operator only.

Jacobian operator will be constructed internally based on the parameter "Jacobian Operator" in the parameter list. See the parameter key "Jacobian Operator" for more details. Defaults to using a NOX::Epetra::MatrixFree object.

References aztecSolverPtr, createJacobianOperator(), getOperatorType(), precPtr, precType, Teuchos::rcp(), reset(), and tmpVectorPtr.

NOX::Epetra::LinearSystemAztecOO::LinearSystemAztecOO ( Teuchos::ParameterList printingParams,
Teuchos::ParameterList linearSolverParams,
const Teuchos::RCP< NOX::Epetra::Interface::Jacobian > &  iJac,
const Teuchos::RCP< Epetra_Operator > &  J,
const Teuchos::RCP< NOX::Epetra::Interface::Preconditioner > &  iPrec,
const Teuchos::RCP< Epetra_Operator > &  M,
const NOX::Epetra::Vector cloneVector,
const Teuchos::RCP< NOX::Epetra::Scaling scalingObject = Teuchos::null 
)

Constructor with user supplied separate objects for the Jacobian (J) and Preconditioner (M). linearSolverParams is the "Linear Solver" sublist of parameter list.

References aztecSolverPtr, getOperatorType(), jacPtr, jacType, precPtr, precType, Teuchos::rcp(), reset(), and tmpVectorPtr.

Member Function Documentation

bool NOX::Epetra::LinearSystemAztecOO::applyJacobian ( const NOX::Epetra::Vector input,
NOX::Epetra::Vector result 
) const
virtual

Applies Jacobian to the given input vector and puts the answer in the result.

Computes

\[ v = J u, \]

where $ J$ is the Jacobian, $ u$ is the input vector, and $ v$ is the result vector. Returns true if successful.

Implements NOX::Epetra::LinearSystem.

References NOX::Epetra::Vector::getEpetraVector().

bool NOX::Epetra::LinearSystemAztecOO::applyJacobianInverse ( Teuchos::ParameterList params,
const NOX::Epetra::Vector input,
NOX::Epetra::Vector result 
)
virtual

Applies the inverse of the Jacobian matrix to the given input vector and puts the answer in result.

Computes

\[ v = J^{-1} u, \]

where $ J$ is the Jacobian, $ u$ is the input vector, and $ v$ is the result vector.

The parameter list contains the linear solver options.

Implements NOX::Epetra::LinearSystem.

References NOX::Utils::Details, Teuchos::ParameterList::get(), NOX::Epetra::Vector::getEpetraVector(), NOX::Epetra::Vector::init(), Teuchos::is_null(), Epetra_RowMatrix::RowMatrixRowMap(), Teuchos::ParameterList::set(), and Teuchos::ParameterList::sublist().

bool NOX::Epetra::LinearSystemAztecOO::applyJacobianTranspose ( const NOX::Epetra::Vector input,
NOX::Epetra::Vector result 
) const
virtual

Applies Jacobian-Transpose to the given input vector and puts the answer in the result.

Computes

\[ v = J^T u, \]

where $ J$ is the Jacobian, $ u$ is the input vector, and $ v$ is the result vector. Returns true if successful.

Implements NOX::Epetra::LinearSystem.

References NOX::Epetra::Vector::getEpetraVector().

bool NOX::Epetra::LinearSystemAztecOO::applyRightPreconditioning ( bool  useTranspose,
Teuchos::ParameterList params,
const NOX::Epetra::Vector input,
NOX::Epetra::Vector result 
) const
virtual

Apply right preconditiong to the given input vector.

Let $ M$ be a right preconditioner for the Jacobian $ J$; in other words, $ M$ is a matrix such that

\[ JM \approx I. \]

Compute

\[ u = M^{-1} v, \]

where $ u$ is the input vector and $ v$ is the result vector.

If useTranspose is true, then the transpose of the preconditioner is applied:

\[ u = {M^{-1}}^T v, \]

The transpose preconditioner is currently only required for Tensor methods.

The parameter list contains the linear solver options.

Implements NOX::Epetra::LinearSystem.

References Teuchos::ParameterList::get(), NOX::Epetra::Vector::getEpetraVector(), NOX::Epetra::Vector::init(), TEUCHOS_TEST_FOR_EXCEPTION, and NOX::Utils::Warning.

bool NOX::Epetra::LinearSystemAztecOO::createPreconditioner ( const NOX::Epetra::Vector x,
Teuchos::ParameterList p,
bool  recomputeGraph 
) const
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.

Implements NOX::Epetra::LinearSystem.

References NOX::Epetra::Vector::getEpetraVector(), Teuchos::is_null(), and NOX::Utils::LinearSolverDetails.

Teuchos::RCP< const Epetra_Operator > NOX::Epetra::LinearSystemAztecOO::getGeneratedPrecOperator ( ) const
virtual

Return preconditioner operator generated and stored in AztecOO.

Note: This should only be called if hasPreconditioner() returns true.

Implements NOX::Epetra::LinearSystem.

NOX::Epetra::LinearSystemAztecOO::OperatorType NOX::Epetra::LinearSystemAztecOO::getOperatorType ( const Epetra_Operator o)
protectedvirtual

Deletes the AztecOO solver object. This is called when the solution vector for the group is changed. The preconditioner is no longer valid so the solver and preconditioner are destroyed by a call to this method.

Returns the type of operator that is passed into the group constructors.

Uses dynamic casting to identify the underlying object type.

Referenced by LinearSystemAztecOO().

NOX::Epetra::LinearSystem::PreconditionerReusePolicyType NOX::Epetra::LinearSystemAztecOO::getPreconditionerPolicy ( bool  advanceReuseCounter = true)
virtual

Evaluates the preconditioner policy at the current state.

NOTE: This can change values between nonlienar iterations. It is not a static value.

Implements NOX::Epetra::LinearSystem.

References NOX::Utils::Details.

bool NOX::Epetra::LinearSystemAztecOO::recomputePreconditioner ( const NOX::Epetra::Vector x,
Teuchos::ParameterList linearSolverParams 
) const
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.

Implements NOX::Epetra::LinearSystem.

References NOX::Epetra::Vector::getEpetraVector(), and NOX::Utils::LinearSolverDetails.

void NOX::Epetra::LinearSystemAztecOO::setAztecOOJacobian ( ) const
protectedvirtual

Sets the epetra Jacobian operator in the AztecOO object.

Makes a call to SetUserMatrix or SetUserOperator to set the Jacobian.

void NOX::Epetra::LinearSystemAztecOO::setAztecOOPreconditioner ( ) const
protectedvirtual

Sets the epetra Preconditioner operator in the AztecOO object.

Makes a call to SetUserOperator. This must be done AFTER the Jacobian matrix/operators is set by setAztecOOJacobian(), otherwise the aztec object may ignore this operation.

References Teuchos::is_null().

void NOX::Epetra::LinearSystemAztecOO::setPrecOperatorForSolve ( const Teuchos::RCP< const Epetra_Operator > &  solvePrecOp)
virtual

Set preconditioner operator for solve.

Note: This should only be called if hasPreconditioner() returns true.

Implements NOX::Epetra::LinearSystem.


The documentation for this class was generated from the following files: