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

Newton direction computation More...

#include <NOX_Direction_Newton.H>

Inheritance diagram for NOX::Direction::Newton:
Inheritance graph
[legend]
Collaboration diagram for NOX::Direction::Newton:
Collaboration graph
[legend]

Public Member Functions

 Newton (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList &params)
 Constructor.
 
virtual ~Newton ()
 Destructor.
 
virtual bool reset (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList &params)
 Reset direction based on possibly new parameters.
 
virtual bool compute (NOX::Abstract::Vector &dir, NOX::Abstract::Group &grp, const NOX::Solver::Generic &solver)
 Compute the direction vector, dir, for a specific method given the current group, grp. More...
 
virtual bool compute (NOX::Abstract::Vector &dir, NOX::Abstract::Group &grp, const NOX::Solver::LineSearchBased &solver)
 Same as compute(NOX::Abstract::Vector&, NOX::Abstract::Group&, const NOX::Solver::Generic&) More...
 
- Public Member Functions inherited from NOX::Direction::Generic
 Generic ()
 Constructor. More...
 
virtual ~Generic ()
 Destructor.
 

Protected Member Functions

virtual bool resetForcingTerm (const NOX::Abstract::Group &soln, const NOX::Abstract::Group &oldSoln, int niter, const NOX::Solver::Generic &solver)
 

Detailed Description

Newton direction computation

Computes the Newton direction by solving the Newton system.

\[ Jd = -F \]

Here $J$ is the n x n Jacobian matrix at the current iterate, $F$ is the n-vector representing the nonlinear function at the current iterate, and $d$ is the n-vector that we are solving for.

If we use an iterative linear solver for the Newton system, then this is called an inexact Newton method. The tolerance used to terminate the linear solve is called the forcing term. The forcing term may be constant, or it may be adjustable. In either case, at iteration $k$ we require,

\[ \frac{\|J_k d_k - (-F_k)\|}{\|F_k\|} \leq \eta_k. \]

Here $\eta_k$ is the forcing term for iteration $k$.

Note
This solution tolerance is to be enforced by the user's implementation of NOX::Abstract::Group::computeNewton; it is passed in as the "Tolerance" in the parameter list for that function.

Adjustable forcing terms were introduced by Eisenstat and Walker (1982); here they are implemented as described in Pernice and Walker (1998). We have two choices for adjustable forcing terms:

Parameters

"Direction":

"Direction"/"Newton":

"Direction"/"Newton"/"Linear Solver":

Note
When using a forcing term, it's critically important the the residual of the original system is used in the comparison. This can be an issue if scaling or left preconditioning is applied to the linear system.

References

Member Function Documentation

bool NOX::Direction::Newton::compute ( NOX::Abstract::Vector dir,
NOX::Abstract::Group grp,
const NOX::Solver::Generic solver 
)
virtual

Compute the direction vector, dir, for a specific method given the current group, grp.

The grp is not const so that we can compute the F vector, the Jacobian matrix, the Newton vector, and so on.

Const access to the solver is used for getting additional information such as the past solution, the iteration number, and so on.

Implements NOX::Direction::Generic.

References NOX::Abstract::Group::computeF(), NOX::Abstract::Group::computeJacobian(), NOX::Abstract::Group::computeNewton(), NOX::Utils::Details, NOX::Abstract::Group::getNewton(), NOX::Solver::Generic::getNumIterations(), NOX::Solver::Generic::getPreviousSolutionGroup(), NOX::Abstract::Group::logLastLinearSolveStats(), NOX::Abstract::Group::Ok, and NOX::Utils::Warning.

bool NOX::Direction::Newton::compute ( NOX::Abstract::Vector dir,
NOX::Abstract::Group grp,
const NOX::Solver::LineSearchBased solver 
)
virtual

Same as compute(NOX::Abstract::Vector&, NOX::Abstract::Group&, const NOX::Solver::Generic&)

Enables direct support for line search based solvers for the purpose of efficiency since the LineSearchBased object has a getStep() function that some directions require.

If it is not redefined in the derived class, it will just call the compute with the NOX::Solver::Generic argument.

Reimplemented from NOX::Direction::Generic.

References NOX::Direction::Generic::compute().

bool NOX::Direction::Newton::resetForcingTerm ( const NOX::Abstract::Group soln,
const NOX::Abstract::Group oldSoln,
int  niter,
const NOX::Solver::Generic solver 
)
protectedvirtual

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