| 
    NOX
    Development
    
   | 
 
Newton-like solver using a trust region. More...
#include <NOX_Solver_TrustRegionBased.H>


Public Member Functions | |
| TrustRegionBased (const Teuchos::RCP< NOX::Abstract::Group > &grp, const Teuchos::RCP< NOX::StatusTest::Generic > &tests, const Teuchos::RCP< Teuchos::ParameterList > ¶ms) | |
| Constructor.  More... | |
| virtual | ~TrustRegionBased () | 
| 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 | StepType { Newton, Cauchy, Dogleg } | 
| Enumerated list for each direction that may be required in the Trust region computation.  More... | |
Protected Member Functions | |
| virtual void | init () | 
| Print out initialization information and calcuation the RHS.  | |
| virtual void | invalid (const std::string ¶m, double value) const | 
| Print and error message and throw and error.  | |
| virtual void | printUpdate () | 
| Prints the current iteration information.  | |
Protected Attributes | |
| Teuchos::RCP< NOX::GlobalData > | globalDataPtr | 
| Pointer to the global data object.  | |
| Teuchos::RCP< NOX::Utils > | utilsPtr | 
| Printing Utils.  | |
| 
Teuchos::RCP < NOX::Abstract::Group >  | solnPtr | 
| Current solution.  | |
| Teuchos::RCP < NOX::Abstract::Group >  | oldSolnPtr | 
| Previous solution pointer.  More... | |
| Teuchos::RCP < NOX::Abstract::Vector >  | newtonVecPtr | 
| Current search direction.pointer.  More... | |
| Teuchos::RCP < NOX::Abstract::Vector >  | cauchyVecPtr | 
| Current search direction.pointer.  More... | |
| Teuchos::RCP < NOX::Abstract::Vector >  | aVecPtr | 
| Extra vector used in computations.  More... | |
| Teuchos::RCP < NOX::Abstract::Vector >  | bVecPtr | 
| Extra vector used in computations.  More... | |
| 
Teuchos::RCP < NOX::StatusTest::Generic >  | testPtr | 
| Stopping test.  | |
| NOX::StatusTest::CheckType | checkType | 
| Type of check to use for status tests. See NOX::StatusTest for more details.  | |
| 
Teuchos::RCP < Teuchos::ParameterList >  | paramsPtr | 
| Input parameters.  | |
| 
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 | initRadius | 
| Initial trust region radius.  | |
| 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.  | |
| NOX::StatusTest::StatusType | status | 
| Status of nonlinear solver.  | |
| StepType | stepType | 
| Type of step to be taken.  | |
| 
Teuchos::RCP < NOX::MeritFunction::Generic >  | meritFuncPtr | 
| Stores a user supplied merit function if supplied in the parameter list.  | |
| bool | useAredPredRatio | 
| If set to true, the minimum improvement ratio condition uses an Ared/Pred approach.  | |
| 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
, 
, 
, 
 is the Jacobian of 
 at 
, and 
 is the trust region radius. 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#245, then shrink the
     trust region: \form#246@_fakenl.
<li> Otherwise if \form#247 and \form#248@_fakenl, then expand the trust region: \form#249@_fakenl.
</ul>
Input Paramters
The following parameters should be specified in the "Trust Region" sublist based to the solver.
"Minimum Trust Region Radius" ( 
) - Minimum allowable trust region radius. Defaults to 1.0e-6.
) - Maximum allowable trust region radius. Defaults to 1.0e+10.
) - Minimum improvement ratio to accept the step. Defaults to 1.0e-4.
) - If the improvement ratio is less than this value, then the trust region is contracted by the amount specified by the "Contraction Factor". Must be larger than "Minimum
  Improvement Ratio". Defaults to 0.1.
) - See above. Defaults to 0.25.
) - If the improvement ratio is greater than this value, then the trust region is contracted by the amount specified by the "Expansion
  Factor". Defaults to 0.75.
) - See above. Defaults to 4.0.
, as described above. The improvement ratio is replaced by an "Ared/Pred" sufficient decrease criteria similar to that used in line search algorithms (see Eisenstat and Walker, SIAM Journal on Optimization V4 no. 2 (1994) pp 393-422):
Output Paramters
A sublist for output parameters called "Output" will be created and contain the following parameters:
"2-Norm or Residual" - Two-norm of final residual
      
  | 
  protected | 
| TrustRegionBased::TrustRegionBased | ( | 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 globalDataPtr, init(), meritFuncPtr, observer, NOX::Solver::parseObserver(), Teuchos::rcp(), Teuchos::ParameterList::sublist(), utilsPtr, and NOX::Solver::validateSolverOptionsSublist().
      
  | 
  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 dx, nIter, status, and NOX::StatusTest::Unconverged.
Referenced by solve().
      
  | 
  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 NOX::Abstract::Group::getNormF(), nIter, observer, paramsPtr, reset(), Teuchos::ParameterList::set(), solnPtr, status, step(), Teuchos::ParameterList::sublist(), and NOX::StatusTest::Unconverged.
      
  | 
  protected | 
Extra vector used in computations.
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
Referenced by step().
      
  | 
  protected | 
Extra vector used in computations.
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
Referenced by step().
      
  | 
  protected | 
Current search direction.pointer.
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
Referenced by step().
      
  | 
  protected | 
Current search direction.pointer.
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
Referenced by step().
      
  | 
  protected | 
Previous solution pointer.
We have both a pointer and a reference because we need to create a DERIVED object and then want to have a reference to it.
Referenced by getPreviousSolutionGroup(), getPreviousSolutionGroupPtr(), and step().
      
  | 
  protected | 
 1.8.5