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::LineSearch::Polynomial Class Reference

A polynomial line search, either quadratic or cubic. More...

#include <NOX_LineSearch_Polynomial.H>

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

Public Member Functions

 Polynomial (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList &params)
 Constructor.
 
 ~Polynomial ()
 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.
 

Protected Types

enum  SufficientDecreaseType { ArmijoGoldstein, AredPred, None }
 Types of sufficient decrease conditions used by checkConvergence() More...
 
enum  InterpolationType { Quadratic, Cubic, Quadratic3 }
 Interpolation types used by compute(). More...
 
enum  RecoveryStepType { Constant, LastComputedStep }
 Type of recovery step to use. More...
 

Protected Member Functions

bool checkConvergence (double newValue, double oldValue, double oldSlope, double step, double eta, int nIters, int nNonlinearIters) const
 Returns true if converged. More...
 
bool updateGrp (NOX::Abstract::Group &newGrp, const NOX::Abstract::Group &oldGrp, const NOX::Abstract::Vector &dir, double step) const
 Updates the newGrp by computing a new x and a new F(x) More...
 
double computeValue (const NOX::Abstract::Group &grp, double phi)
 Compute the value used to determine sufficient decrease. More...
 
void printOpeningRemarks () const
 Used to print opening remarks for each call to compute().
 
void printBadSlopeWarning (double slope) const
 Prints a warning message saying that the slope is negative.
 

Protected Attributes

SufficientDecreaseType suffDecrCond
 Choice for sufficient decrease condition; uses "Sufficient Decrease Condition" parameter.
 
InterpolationType interpolationType
 Choice of interpolation type; uses "Interpolation Type" parameter.
 
RecoveryStepType recoveryStepType
 Choice of the recovery step type; uses "Recovery Step Type" parameter.
 
double minStep
 Minimum step length (i.e., when we give up); uses "Mimimum Step" parameter.
 
double defaultStep
 Default step; uses "Default Step" parameter.
 
double recoveryStep
 Default step for linesearch failure; uses "Recovery Step" parameter.
 
int maxIters
 Maximum iterations; uses "Max Iters" parameter.
 
double alpha
 The $ \alpha $ for the Armijo-Goldstein condition, or the $ \alpha $ for the Ared/Pred condition; see checkConvergence(). Uses "Alpha Factor" parameter.
 
double minBoundFactor
 Factor that limits the minimum size of the new step based on the previous step; uses "Min Bounds Factor" parameter.
 
double maxBoundFactor
 Factor that limits the maximum size of the new step based on the previous step; uses "Max Bounds Factor" parameter.
 
bool doForceInterpolation
 True is we should force at least one interpolation step; uses "Force Interpolation" parameter.
 
int maxIncreaseIter
 No increases are allowed if the number of nonlinear iterations is greater than this value; uses "Maximum Iterations for Increase".
 
bool doAllowIncrease
 True if we sometimes allow an increase(!) in the decrease measure, i.e., maxIncreaseIter > 0.
 
double maxRelativeIncrease
 Maximum allowable relative increase for decrease meausre (if allowIncrease is true); uses "Allowed Relative Increase" parameter.
 
bool useCounter
 True if we should use counter and output the results; uses "Use Counters" parameter.
 
Teuchos::RCP< NOX::GlobalDataglobalDataPtr
 Global data pointer. Keep this so the parameter list remains valid.
 
Teuchos::ParameterListparamsPtr
 Pointer to the input parameter list. More...
 
NOX::LineSearch::Utils::Printing print
 Common line search printing utilities.
 
NOX::LineSearchCounterscounter
 Common common counters for line searches.
 
NOX::LineSearch::Utils::Slope slopeUtil
 Common slope calculations for line searches.
 
Teuchos::RCP
< NOX::MeritFunction::Generic
meritFuncPtr
 Pointer to a user supplied merit function.
 

Detailed Description

A polynomial line search, either quadratic or cubic.

This line search can be called via NOX::LineSearch::Manager.

The goal of the line search is to find a step $\lambda$ for the calculation $x_{\rm new} = x_{\rm old} + \lambda d$, given $x_{\rm old}$ and $ d $. To accomplish this goal, we minimize a merit function $ \phi(\lambda) $ that measures the "goodness" of the step $\lambda$. The standard merit function is

\[ \phi(\lambda) \equiv \frac{1}{2}||F (x_{\rm old} + \lambda s)||^2, \]

but a user defined merit function can be used instead (see computePhi() for details). Our first attempt is always to try a step $ \lambda_0 $, and then check the stopping criteria. The value of $ \lambda_0 $ is the specified by the "Default Step" parameter. If the first try doesn't work, then we successively minimize polynomial approximations, $ p_k(\lambda) \approx \phi(\lambda) $.

Stopping Criteria

The inner iterations continue until:

Polynomial Models of the Merit Function

We compute $ p_k(\lambda) $ by interpolating $f$. For the quadratic fit, we interpolate $ \phi(0) $, $ \phi'(0) $, and $ \phi(\lambda_{k-1}) $ where $ \lambda_{k-1} $ is the $ k-1 $ approximation to the step. For the cubit fit, we additionally include $\phi(\lambda{k-2})$.

The steps are calculated iteratively as follows, depending on the choice of "Interpolation Type".

For "Quadratic" and "Cubic", see computeSlope() for details on how $ \phi'(0) $ is calculated.

Bounds on the step length

We do not allow the step to grow or shrink too quickly by enforcing the following bounds:

\[ \gamma_{min} \; \lambda_{k-1} \leq \lambda_k \le \gamma_{max} \; \lambda_{k-1} \]

Here $ \gamma_{min} $ and $ \gamma_{max} $ are defined by parameters "Min Bounds Factor" and "Max Bounds Factor".

Input Parameters

"Line Search":

"Line Search"/"Polynomial":

Output Parameters

If the "Use Counters" parameter is set to true, then a sublist for output parameters called "Output" will be created in the parameter list used to instantiate or reset the class. Valid output parameters are:

References

This line search is based on materials in the following:

Author
Russ Hooper, Roger Pawlowski, Tammy Kolda

Member Enumeration Documentation

Interpolation types used by compute().

Enumerator
Quadratic 

Use quadratic interpolation throughout.

Cubic 

Use quadratic interpolation in the first inner iteration and cubic interpolation otherwise.

Quadratic3 

Use a 3-point quadratic line search. (The second step is 0.5 times the default step.)

Type of recovery step to use.

Enumerator
Constant 

Use a constant value.

LastComputedStep 

Use the last value computed in the line search algorithm.

Types of sufficient decrease conditions used by checkConvergence()

Enumerator
ArmijoGoldstein 

Armijo-Goldstein conditions.

AredPred 

Ared/Pred condition.

None 

No condition.

Member Function Documentation

bool NOX::LineSearch::Polynomial::checkConvergence ( double  newValue,
double  oldValue,
double  oldSlope,
double  step,
double  eta,
int  nIters,
int  nNonlinearIters 
) const
protected

Returns true if converged.

We go through the following checks for convergence.

  1. If the "Force Interpolation" parameter is true and this is the first inner iteration (i.e., nIters is 1), then we return false.

  2. Next, it checks to see if the "Relative Increase" condition is satisfied. If so, then we return true. The "Relative Increase" condition is satisfied if both of the following two conditions hold.

    • The number of nonlinear iterations is less than or equal to the value specified in the "Maximum Iteration for Increase" parameter
    • The ratio of newValue to oldValue is less than the value specified in "Allowed Relative Increase".

  3. Last, we check the sufficient decrease condition. We return true if it's satisfied, and false otherwise. The condition depends on the value of "Sufficient Decrease Condition", as follows.

    • "Armijo-Goldstein": Return true if

      \[ \phi(\lambda) \leq \phi(0) + \alpha \; \lambda \; \phi'(0) \]

    • "Ared/Pred": Return true if

      \[ \| F(x_{\rm old} + \lambda d )\| \leq \| F(x_{\rm old}) \| (1 - \alpha (1 - \eta)) \]

    • "None": Always return true.

    For the first two cases, $\alpha$ is specified by the parameter "Alpha Factor".

Parameters
newValueDepends on the "Sufficient Decrease Condition" parameter.
  • "Armijo-Goldstein" - $ \phi(\lambda) $
  • "Ared/Pred" - $ \| F(x_{\rm old} + \lambda d )\| $
  • "None" - N/A
oldValueDepends on the "Sufficient Decrease Condition" parameter.
  • "Armijo-Goldstein" - $ \phi(0) $
  • "Ared/Pred" - $ \| F(x_{\rm old} )\| $
  • "None" - N/A
oldSlopeOnly applicable to "Armijo-Goldstein", in which case it should be $\phi'(0)$.
stepThe current step, $\lambda$.
etaOnly applicable to "Ared/Pred", in which case it should be the value of $\eta$ for last forcing term calculation in NOX::Direction::Newton
nItersNumber of inner iterations.
nNonlinearItersNumber of nonlinear iterations.
Note
The norm used for "Ared/Pred" can be specified by the user by using the "User Defined Norm" parameter; this parameter is any object derived from NOX::Parameter::UserNorm.

References NOX::StatusTest::FiniteValue::finiteNumberTest().

bool NOX::LineSearch::Polynomial::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 Teuchos::ParameterList::get(), NOX::Solver::Generic::getList(), NOX::Solver::Generic::getNumIterations(), and NOX::Solver::Generic::getPreviousSolutionGroup().

double NOX::LineSearch::Polynomial::computeValue ( const NOX::Abstract::Group grp,
double  phi 
)
protected

Compute the value used to determine sufficient decrease.

If the "Sufficient Decrease Condition" is set to "Ared/Pred", then we do the following. If a "User Defined Norm" parameter is specified, then the NOX::Parameter::UserNorm::norm function on the user-supplied merit function is used. If the user does not provide a norm, we return $ \|F(x)\| $.

If the "Sufficient Decrease Condition" is not set to "Ared/Pred", then we simply return phi.

Parameters
phi- Should be equal to computePhi(grp).

References NOX::Abstract::Group::getNormF().

bool NOX::LineSearch::Polynomial::updateGrp ( NOX::Abstract::Group newGrp,
const NOX::Abstract::Group oldGrp,
const NOX::Abstract::Vector dir,
double  step 
) const
protected

Updates the newGrp by computing a new x and a new F(x)

Let

  • $x_{\rm new}$ denote the new solution to be calculated (corresponding to newGrp)
  • $x_{\rm old}$ denote the previous solution (i.e., the result of oldGrp.getX())
  • $d$ denotes the search direction (dir).
  • $\lambda$ denote the step (step),

We compute

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

and $ F(x_{\rm new})$. The results are stored in newGrp.

References NOX::Abstract::Group::computeF(), NOX::Abstract::Group::computeX(), and NOX::Abstract::Group::Ok.

Member Data Documentation

Teuchos::ParameterList* NOX::LineSearch::Polynomial::paramsPtr
protected

Pointer to the input parameter list.

We need this to later create an "Output" sublist to store output parameters from counter.


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