NOX
Development
|
Newton-like solver using a trust region. More...
#include <NOX_Solver_InexactTrustRegionBased.H>
Public Member Functions | |
InexactTrustRegionBased (const Teuchos::RCP< NOX::Abstract::Group > &grp, const Teuchos::RCP< NOX::StatusTest::Generic > &tests, const Teuchos::RCP< Teuchos::ParameterList > ¶ms) | |
Constructor. More... | |
virtual | ~InexactTrustRegionBased () |
Destructor. | |
virtual void | reset (const NOX::Abstract::Vector &initialGuess, const Teuchos::RCP< NOX::StatusTest::Generic > &tests) |
Resets the solver, sets a new status test, and sets a new initial guess. | |
virtual void | reset (const NOX::Abstract::Vector &initialGuess) |
Resets the solver and sets a new initial guess. | |
virtual void | reset () |
Resets the solver for another solve. This resets the counters and status only. Uses the final solution from the last solve as the initial guess for the next solve. More... | |
virtual NOX::StatusTest::StatusType | step () |
Do one nonlinear step in the iteration sequence and return status. | |
virtual NOX::StatusTest::StatusType | solve () |
Solve the nonlinear problem and return final status. More... | |
virtual const NOX::Abstract::Group & | getSolutionGroup () const |
Return a reference to the current solution group. | |
virtual const NOX::Abstract::Group & | getPreviousSolutionGroup () const |
Return a reference to the previous solution group. | |
virtual NOX::StatusTest::StatusType | getStatus () const |
Returns the current status of the solver. | |
virtual int | getNumIterations () const |
Get number of iterations. | |
virtual const Teuchos::ParameterList & | getList () const |
Return a reference to the solver parameters. | |
virtual Teuchos::RCP< const NOX::Abstract::Group > | getSolutionGroupPtr () const |
Return a RCP to the solution group. | |
virtual Teuchos::RCP< const NOX::Abstract::Group > | getPreviousSolutionGroupPtr () const |
Return a RCP to the previous solution group. | |
virtual Teuchos::RCP< const Teuchos::ParameterList > | getListPtr () const |
Return a RCP to the solver parameters. | |
virtual Teuchos::RCP< const NOX::SolverStats > | getSolverStatistics () const |
Return a RCP to the solver statistics. | |
Public Member Functions inherited from NOX::Solver::Generic | |
Generic () | |
Constructor (does nothing) | |
virtual | ~Generic () |
Destructor (does nothing) | |
Protected Types | |
enum | TrustRegionType { Standard, Inexact } |
Type of Trust Region algorithm to use. More... | |
enum | InnerIterationReturnType { Converged, Unconverged, Failed } |
Return types for inner iteration status test. More... | |
enum | StepType { Newton, Cauchy, Dogleg } |
Enumerated list for each direction that may be required in the Trust region computation. More... | |
Protected Member Functions | |
virtual NOX::StatusTest::StatusType | iterateStandard () |
"Standard" trust region implementation | |
virtual NOX::StatusTest::StatusType | iterateInexact () |
"Inexact Trust Region" | |
virtual void | init () |
Print out initialization information and calcuation the RHS. | |
virtual void | printUpdate () |
Prints the current iteration information. | |
virtual void | invalid (const std::string ¶m, double value) const |
Print an error message and throw an error during parameter reads. | |
virtual void | throwError (const std::string &method, const std::string &mesage) const |
Print an error message and throw an error. | |
NOX::StatusTest::StatusType | checkStep (const NOX::Abstract::Vector &step, double &radius) |
Check to see if the current step is acceptable. If not, it reduces the trust region radius accordingly. | |
virtual double | computeNorm (const NOX::Abstract::Vector &v) |
Computes the norm of a given vector. More... | |
Protected Attributes | |
TrustRegionType | method |
Type of trust region algorithm to use. | |
Teuchos::RCP< NOX::GlobalData > | globalDataPtr |
Pointer to the global data object. | |
Teuchos::RCP< NOX::Utils > | utils |
Utils. | |
InnerIterationReturnType | innerIterationStatus |
Current status of the trust region inner iteration. | |
Teuchos::RCP < NOX::Abstract::Group > | solnPtr |
Current solution. | |
Teuchos::RCP < NOX::Abstract::Group > | oldSolnPtr |
Previous solution pointer. | |
Teuchos::RCP < NOX::Abstract::Vector > | newtonVecPtr |
Current newton direction pointer. | |
Teuchos::RCP < NOX::Abstract::Vector > | cauchyVecPtr |
Current cauchy direction pointer. | |
Teuchos::RCP < NOX::Abstract::Vector > | rCauchyVecPtr |
Extra vector used in computations. | |
Teuchos::RCP < NOX::Abstract::Vector > | residualVecPtr |
Extra vector used in computations. | |
Teuchos::RCP < NOX::Abstract::Vector > | aVecPtr |
Extra vector used in computations. | |
Teuchos::RCP < NOX::Abstract::Vector > | bVecPtr |
Extra vector used in computations. | |
Teuchos::RCP < NOX::StatusTest::Generic > | testPtr |
Stopping test. | |
Teuchos::RCP < Teuchos::ParameterList > | paramsPtr |
Input parameters. | |
NOX::Direction::Utils::InexactNewton | inNewtonUtils |
Inexact Newton utitilities. | |
Teuchos::RCP < NOX::Direction::Generic > | newtonPtr |
Newton Search Direction. | |
Teuchos::RCP < NOX::Direction::Generic > | cauchyPtr |
Cauchy Search Direction. | |
double | radius |
Radius of the trust region. | |
double | minRatio |
Minimum improvement ratio to accept step. | |
double | minRadius |
Minimum trust region radius. | |
double | maxRadius |
Maximum trust region radius. | |
double | contractTriggerRatio |
ratio < alpha triggers contraction | |
double | expandTriggerRatio |
ratio > beta triggers expansion | |
double | expandFactor |
Expansion factor. | |
double | contractFactor |
Constraction factor. | |
double | recoveryStep |
double | newF |
Value of at current solution. | |
double | oldF |
Value of at previous solution. | |
double | dx |
norm(xnew - xold) | |
int | nIter |
Number of nonlinear iterations. | |
double | eta |
Current linear solve tolerance (inexact only). | |
double | eta_last |
Linear solve tolerance used in last iteration (inexact only). | |
NOX::StatusTest::StatusType | status |
Status of nonlinear solver. | |
NOX::StatusTest::CheckType | checkType |
Type of check to use for status tests. See NOX::StatusTest for more details. | |
StepType | stepType |
Type of step to be taken. | |
Teuchos::RCP < NOX::MeritFunction::Generic > | meritFuncPtr |
Stores merit function supplied by global data. | |
bool | useCauchyInNewtonDirection |
If set to true, the initial guess for the Newton direction computation will use the Cauchy direction as the initial guess. | |
bool | writeOutputParamsToList |
If set to true, statistics/counters will be output to the output list. | |
bool | useCounters |
If set to true, counters will be stored by the solver. | |
NOX::SolverStats::TrustRegionStats * | counters |
Counters for the algorithm. | |
bool | useAredPredRatio |
If set to true, the minimum improvement ratio condition uses an Ared/Pred approach. | |
bool | useDoglegMinimization |
If set to true, the parameter is minimized over the dogleg line segments instead of being computed at the trust regioin radius. | |
Teuchos::RCP< NOX::Observer > | observer |
Pointer to a user defined NOX::Observer object. | |
Newton-like solver using a trust region.
Our goal is to solve: where . Alternatively, we might say that we wish to solve
The trust region subproblem (TRSP) at iteration is given by
where
The "improvement ratio" for a given step is defined as
An iteration consists of the following steps.
Compute Newton-like direction:
Compute Cauchy-like direction:
If this is the first iteration, initialize as follows: If , then ; else, .
Initialize
While and , do the following.
Compute the direction as follows:
If , then take a Newton step by setting
Otherwise if , then take a Cauchy step by setting
Otherwise, take a Dog Leg step by setting where with .
Set and calculate
If , then Otherwise
</ul>
Update the solution:
Update trust region:
If and , then shrink the trust region to the size of the Newton step: .
<li> Otherwise if \form#238, then shrink the trust region: \form#239@_fakenl. <li> Otherwise if \form#240 and \form#241@_fakenl, then expand the trust region: \form#242@_fakenl. </ul>
Input Paramters
The following parameters should be specified in the "Trust Region" sublist based to the solver.
Output Paramters
A sublist called "Output" will be created at the top level of the parameter list and contain the following general solver parameters:
A sublist called "Output" will be created in the "Trust Region" sublist and contain the following trust region specific output parameters:
"Dogleg Steps: Average Fraction Between Cauchy and Newton Direction" - Average value of the fraction a dogleg step took between the Cauchy and Newton directions. This is the variable in the standard dogleg algorithm and the parameter in the inexact dogleg algorithm. A value of 0.0 is a full step in the Cauchy direction and a value of 1.0 is a full step in the Newton direction.
|
protected |
|
protected |
NOX::Solver::InexactTrustRegionBased::InexactTrustRegionBased | ( | const Teuchos::RCP< NOX::Abstract::Group > & | grp, |
const Teuchos::RCP< NOX::StatusTest::Generic > & | tests, | ||
const Teuchos::RCP< Teuchos::ParameterList > & | params | ||
) |
Constructor.
See reset() for description.
References counters, globalDataPtr, init(), inNewtonUtils, meritFuncPtr, observer, paramsPtr, NOX::Solver::parseObserver(), Teuchos::rcp(), NOX::Direction::Utils::InexactNewton::reset(), Teuchos::ParameterList::sublist(), utils, and NOX::Solver::validateSolverOptionsSublist().
|
protectedvirtual |
Computes the norm of a given vector.
Defaults to the L-2 norm but could use a user defined norm also.
References NOX::Abstract::Vector::norm().
|
virtual |
Resets the solver for another solve. This resets the counters and status only. Uses the final solution from the last solve as the initial guess for the next solve.
NOTE: All NOX solvers will call reset() automatically at teh beginning of the solve() method. We add the reset() method to the solver interface for the application to call in case the application needs to reset counters and status manually before the next call to solve() is made.
Implements NOX::Solver::Generic.
References NOX::Utils::fill(), NOX::Utils::Parameters, and NOX::StatusTest::Unconverged.
|
virtual |
Solve the nonlinear problem and return final status.
By "solve", we call iterate() until the NOX::StatusTest value is either NOX::StatusTest::Converged or NOX::StatusTest::Failed.
Implements NOX::Solver::Generic.
References Teuchos::ParameterList::set(), Teuchos::ParameterList::sublist(), and NOX::StatusTest::Unconverged.
|
protected |
Take a step of this length in the Newton direction if the trust-region search fails