NOX
Development
|
A polynomial line search, either quadratic or cubic. More...
#include <NOX_LineSearch_Polynomial.H>
Public Member Functions | |
Polynomial (const Teuchos::RCP< NOX::GlobalData > &gd, Teuchos::ParameterList ¶ms) | |
Constructor. | |
~Polynomial () | |
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. | |
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 for the Armijo-Goldstein condition, or the 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::GlobalData > | globalDataPtr |
Global data pointer. Keep this so the parameter list remains valid. | |
Teuchos::ParameterList * | paramsPtr |
Pointer to the input parameter list. More... | |
NOX::LineSearch::Utils::Printing | |
Common line search printing utilities. | |
NOX::LineSearchCounters * | counter |
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. | |
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 for the calculation , given and . To accomplish this goal, we minimize a merit function that measures the "goodness" of the step . The standard merit function is
but a user defined merit function can be used instead (see computePhi() for details). Our first attempt is always to try a step , and then check the stopping criteria. The value of is the specified by the "Default Step" parameter. If the first try doesn't work, then we successively minimize polynomial approximations, .
Stopping Criteria
The inner iterations continue until:
The sufficient decrease condition is met; see checkConvergence() for details.
The maximum iterations are reached; see parameter "Max Iters". This is considered a failure and the recovery step is used; see parameter "Recovery Step".
The minimum step length is reached; see parameter "Minimum Step". This is considered a line search failure and the recovery step is used; see parameter "Recovery Step".
Polynomial Models of the Merit Function
We compute by interpolating . For the quadratic fit, we interpolate , , and where is the approximation to the step. For the cubit fit, we additionally include .
The steps are calculated iteratively as follows, depending on the choice of "Interpolation Type".
"Quadratic" - We construct a quadratic model of , and solve for to get
"Cubic" - We construct a cubic model of , and solve for . We use the quadratic model to solve for ; otherwise,
where
"Quadratic3" - We construct a quadratic model of using , , and . No derivative information for is used. We let , and otherwise
For "Quadratic" and "Cubic", see computeSlope() for details on how is calculated.
Bounds on the step length
We do not allow the step to grow or shrink too quickly by enforcing the following bounds:
Here and are defined by parameters "Min Bounds Factor" and "Max Bounds Factor".
Input Parameters
"Line Search":
"Line Search"/"Polynomial":
"Default Step" - Starting step length, i.e., . Defaults to 1.0.
"Max Iters" - Maximum number of line search iterations. The search fails if the number of iterations exceeds this value. Defaults to 100.
"Minimum Step" - Minimum acceptable step length. The search fails if the computed is less than this value. Defaults to 1.0e-12.
"Recovery Step Type" - Determines the step size to take when the line search fails. Choices are:
"Recovery Step" - The value of the step to take when the line search fails. Only used if the "Recovery Step Type" is set to "Constant". Defaults to value for "Default Step".
"Interpolation Type" - Type of interpolation that should be used. Choices are
"Min Bounds Factor" - Choice for , i.e., the factor that limits the minimum size of the new step based on the previous step. Defaults to 0.1.
"Max Bounds Factor" - Choice for , i.e., the factor that limits the maximum size of the new step based on the previous step. Defaults to 0.5.
"Sufficient Decrease Condition" - See checkConvergence() for details. Choices are:
"None"
"Alpha Factor" - Parameter choice for sufficient decrease condition. See checkConvergence() for details. Defaults to 1.0e-4.
"Force Interpolation" (boolean) - Set to true if at least one interpolation step should be used. The default is false which means that the line search will stop if the default step length satisfies the convergence criteria. Defaults to false.
"Use Counters" (boolean) - Set to true if we should use counters and then output the result to the paramter list as described in Output Parameters. Defaults to true.
"Maximum Iteration for Increase" - Maximum index of the nonlinear iteration for which we allow a relative increase. See checkConvergence() for further details. Defaults to 0 (zero).
"Allowed Relative Increase" - See checkConvergence() for details. Defaults to 100.
"User Defined Merit Function" - The user can redefine the merit function used; see computePhi() and NOX::Parameter::MeritFunction for details.
"User Defined Norm" - The user can redefine the norm that is used in the Ared/Pred sufficient decrease condition; see computeValue() and NOX::Parameter::UserNorm for details.
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:
"Total Number of Line Search Calls" - Total number of calls to the compute() method of this line search.
"Total Number of Non-trivial Line Searches" - Total number of steps that could not directly take a full step and meet the required "Sufficient Decrease Condition" (i.e., the line search had to do at least one interpolation).
"Total Number of Failed Line Searches" - Total number of line searches that failed and used a recovery step.
"Total Number of Line Search Inner Iterations" - Total number of inner iterations of all calls to compute().
References
This line search is based on materials in the following:
|
protected |
Interpolation types used by compute().
|
protected |
|
protected |
Types of sufficient decrease conditions used by checkConvergence()
Enumerator | |
---|---|
ArmijoGoldstein |
Armijo-Goldstein conditions. |
AredPred |
Ared/Pred condition. |
None |
No condition. |
|
protected |
Returns true if converged.
We go through the following checks for convergence.
If the "Force Interpolation" parameter is true and this is the first inner iteration (i.e., nIters is 1), then we return false.
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.
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
"Ared/Pred": Return true if
For the first two cases, is specified by the parameter "Alpha Factor".
newValue | Depends on the "Sufficient Decrease Condition" parameter.
|
oldValue | Depends on the "Sufficient Decrease Condition" parameter.
|
oldSlope | Only applicable to "Armijo-Goldstein", in which case it should be . |
step | The current step, . |
eta | Only applicable to "Ared/Pred", in which case it should be the value of for last forcing term calculation in NOX::Direction::Newton |
nIters | Number of inner iterations. |
nNonlinearIters | Number of nonlinear iterations. |
References NOX::StatusTest::FiniteValue::finiteNumberTest().
|
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 Teuchos::ParameterList::get(), NOX::Solver::Generic::getList(), NOX::Solver::Generic::getNumIterations(), and NOX::Solver::Generic::getPreviousSolutionGroup().
|
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 .
If the "Sufficient Decrease Condition" is not set to "Ared/Pred", then we simply return phi.
phi | - Should be equal to computePhi(grp). |
References NOX::Abstract::Group::getNormF().
|
protected |
Updates the newGrp by computing a new x and a new F(x)
Let
newGrp
)dir
).step
),We compute
and . The results are stored in newGrp
.
References NOX::Abstract::Group::computeF(), NOX::Abstract::Group::computeX(), and NOX::Abstract::Group::Ok.
|
protected |
Pointer to the input parameter list.
We need this to later create an "Output" sublist to store output parameters from counter.