NOX
Development
|
Use NonlinearCG linesearch. More...
#include <NOX_LineSearch_NonlinearCG.H>
Public Member Functions | |
NonlinearCG (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList ¶ms) | |
Constructor. | |
~NonlinearCG () | |
Destructor. | |
bool | reset (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList ¶ms) |
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. | |
Use NonlinearCG linesearch.
This is a simple linesearch intended to be used with NOX::Direction::NonlinearCG, which provides search direction , in computing an update to the current solution vector . It is designed to compute a step length 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, is computed from a single evaluation of,
where is the n x n Jacobian matrix. Explicit construction of is avoided by performing the product using a directional derivative (cf NOX::Epetra::MatrixFree):
where .
Derivation / Theory:
This linesearch is derived by attempting to achieve in a single step, the following minimization:
where is a merit function chosen (but never explicitly given) so that an equivalence to Linear CG holds, ie . The minimization above can now be cast as an equation:
An approximate solution to this equation can be obtained from a second-order expansion of
,
from which it immediately follows
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 .
|
virtual |
Perform a line search.
On input:
grp | The initial solution vector, . |
dir | A vector of directions to be used in the line search, . |
s | The nonlinear solver. |
On output:
step | The distance the direction was scaled, . |
grp | The grp is updated with a new solution, , resulting from the linesearch. Normally, for a single direction line search, this is computed as: |
Ideally, (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().