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

Use NonlinearCG linesearch. More...

#include <NOX_LineSearch_NonlinearCG.H>

Inheritance diagram for NOX::LineSearch::NonlinearCG:
Inheritance graph
[legend]
Collaboration diagram for NOX::LineSearch::NonlinearCG:
Collaboration graph
[legend]

Public Member Functions

 NonlinearCG (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList &params)
 Constructor.
 
 ~NonlinearCG ()
 Destructor.
 
bool reset (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList &params)
 
bool compute (NOX::Abstract::Group &newgrp, double &step, const NOX::Abstract::Vector &dir, const NOX::Solver::Generic &s)
 Perform a line search. More...
 
- Public Member Functions inherited from NOX::LineSearch::Generic
 Generic ()
 Default constructor.
 
virtual ~Generic ()
 Destructor.
 

Detailed Description

Use NonlinearCG linesearch.

This is a simple linesearch intended to be used with NOX::Direction::NonlinearCG, which provides search direction $ d $, in computing an update to the current solution vector $ x_{new} = x_{old} + \lambda d $. It is designed to compute a step length $ \lambda $ consistent with the exact linesearch of Linear CG for linear problems, and it avoids use of matrices by employing a directional derivative (details below). The step length, $ \lambda $ is computed from a single evaluation of,

\[ \lambda = - \frac{F(x_{old})^T d}{d^T J(x_{old})d} \]

where $ J $ is the n x n Jacobian matrix. Explicit construction of $ J $ is avoided by performing the product $ Jd $ using a directional derivative (cf NOX::Epetra::MatrixFree):

\[ J(x_{old})d \approx \frac{F(x_{old} + \delta d) - F(x_{old})}{\delta} \]

where $ \delta = 10^{-6} (10^{-6} + ||x_{old}|| / ||d||) $ .

Derivation / Theory:

This linesearch is derived by attempting to achieve in a single step, the following minimization:

\[ \min_\lambda \phi(\lambda)\equiv\phi (x_{old}+ \lambda d) \]

where $ \phi $ is a merit function chosen (but never explicitly given) so that an equivalence to Linear CG holds, ie $ \nabla\phi(x) \leftrightarrow F(x) $. The minimization above can now be cast as an equation:

\[ \phi ' (\lambda) = \nabla\phi (x_{old}+ \lambda d)^T d = F(x_{old}+ \lambda d)^T d = 0~~. \]

An approximate solution to this equation can be obtained from a second-order expansion of

\[ \phi(\lambda) \]

,

\[ \phi(\lambda)\approx\phi (0) + \phi ' (0)\lambda + \phi '' (0) \frac{\lambda^2}{2} \]

from which it immediately follows

\[ \%lambda_{min} \approx - \frac{\phi ' (0)}{\phi '' (0)} = - \frac{F(x_{old})^T d}{d^T J(x_{old})d} \]

Input Parameters

The NonlinearCG linesearch is selected using:

"Line Search"

"Method" = "%NonlinearCG" [required]

Currently, no adjustable parameters exist for this linesarch.

References

linesearch is adapted from ideas presented in Section 14.2 of:

Jonathan Richard Shewchuk, "An Introduction to the Conjugate Gradient Method Without the Agonizing Pain," 1994. Though presented within the context of nonlinear optimization, the connection to solving nonlinear equation systems is made via the equivalence $ f'(x) \leftrightarrow F(x) $.

Author
Russ Hooper, Org. 9233, Sandia National Labs

Member Function Documentation

bool NOX::LineSearch::NonlinearCG::compute ( NOX::Abstract::Group grp,
double &  step,
const NOX::Abstract::Vector dir,
const NOX::Solver::Generic s 
)
virtual

Perform a line search.

On input:

Parameters
grpThe initial solution vector, $x_{\rm old}$.
dirA vector of directions to be used in the line search, $d$.
sThe nonlinear solver.

On output:

Parameters
stepThe distance the direction was scaled, $ \lambda $.
grpThe grp is updated with a new solution, $ x_{\rm new} $, resulting from the linesearch. Normally, for a single direction line search, this is computed as:

\[ x_{\rm new} = x_{\rm old} + \lambda d. \]

Ideally, $ \|F(x_{\rm new})\| < \|F(x_{\rm old})\| $ (e.g the final direction is a descent direction).

Note that the dir object is a std::vector. For typical line searches as described in the above equation, this vector is of size one. We have used a std::vector to allow for special cases of multi-directional line searches such as the Bader/Schnabel curvillinear line search.

Return value is true for a successful line search computation.

Implements NOX::LineSearch::Generic.

References NOX::Abstract::Group::computeF(), NOX::Abstract::Group::computeX(), NOX::Utils::fill(), NOX::Abstract::Group::getF(), NOX::Solver::Generic::getPreviousSolutionGroup(), NOX::Utils::InnerIteration, and NOX::Abstract::Vector::innerProduct().


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