OPTPP::OptNewtonLike Class Reference

#include <OptNewtonLike.h>

Inheritance diagram for OPTPP::OptNewtonLike:

OPTPP::OptimizeClass OPTPP::OptNewton1Deriv OPTPP::OptNewton2Deriv OPTPP::OptFDNewton OPTPP::OptQNewton OPTPP::OptNewton List of all members.

Public Member Functions

 OptNewtonLike ()
 OptNewtonLike (int n)
 OptNewtonLike (int n, UPDATEFCN u)
 OptNewtonLike (int n, TOLS t)
virtual ~OptNewtonLike ()
virtual void acceptStep (int k, int step_type)
virtual NEWMAT::ColumnVector computeSearch (NEWMAT::SymmetricMatrix &H)
 Compute Newton direction.
virtual void updateModel (int k, int ndim, NEWMAT::ColumnVector x)
virtual int checkConvg ()
 Check to see if algorithm satisfies the convergence criterion.
virtual int checkDeriv ()
 Compare the analytic gradient with the finite difference gradient.
virtual int computeStep (NEWMAT::ColumnVector sk)
 Compute the step length along the Newton direction.
virtual void initOpt ()
 Initialize algorithmic parameters.
virtual void initHessian ()
 Compute the Hessian or its approximation at the initial point.
virtual double initTrustRegionSize () const
 the trustregion or trustpds globalization strategies are selected
virtual void optimize ()
 Invoke Newton's method on an unconstrained problem.
virtual void readOptInput ()
 Read user-specified input options from a file.
virtual void reset ()
 Reset parameters.
int checkAnalyticFDGrad ()
 Compare the analytic gradient with the finite difference gradient.
int getFevals () const
int getGevals () const
real getTRSize () const
void setTRSize (real delta)
 Set trust-region radius.
real getGradMult () const
void setGradMult (real tau)
 Set gradient multiplier which is used to compute trust-region radius.
int getSearchSize () const
void setSearchSize (int sss)
 Set number of points in search scheme for trust-pds search strategy.
bool getWarmStart () const
void UseWarmStart (NEWMAT::SymmetricMatrix &H)
void setSearchStrategy (SearchStrategy s)
SearchStrategy getSearchStrategy () const
 Set globalization strategy for optimization algorithms.
void setDerivOption (DerivOption d)
DerivOption getDerivOption () const
 Set the type of finite difference routine.
NEWMAT::SymmetricMatrix getHessian () const
void setHessian (NEWMAT::SymmetricMatrix &H)
 Store the current Hessian matrix.
virtual NEWMAT::SymmetricMatrix updateH (NEWMAT::SymmetricMatrix &H, int k)=0
 Compute the Hessian of the objective function or its approximation at the current point.
void printStatus (char *)
 Print status of unconstrained Newton's method.

Protected Member Functions

virtual NLP1nlprob () const =0
 Returns pointer to an NLP1 object.
void defaultAcceptStep (int, int)
NEWMAT::ColumnVector defaultComputeSearch (NEWMAT::SymmetricMatrix &)

Protected Attributes

NEWMAT::ColumnVector gprev
 Gradient at the prev. iteration.
NEWMAT::SymmetricMatrix Hessian
 Current Hessian.
int grad_evals
 Number of gradient evaluations.
SearchStrategy strategy
 User-specified globalization Strategy.
DerivOption finitediff
 User-specified derivative option.
real TR_size
 Trust region radius.
real gradMult
 Gradient multiplier to compute TR_size.
int searchSize
 Search pattern size for TRPDS.
bool WarmStart

Friends

int trustregion (NLP1 *, ostream *, NEWMAT::SymmetricMatrix &, NEWMAT::ColumnVector &, NEWMAT::ColumnVector &, real &, real &, real stpmax, real stpmin)
int trustpds (NLP1 *, ostream *, NEWMAT::SymmetricMatrix &, NEWMAT::ColumnVector &, NEWMAT::ColumnVector &, real &, real &, real stpmax, real stpmin, int)

Detailed Description

OptNewtonLike is the base class for Newton Methods. OptNewtonLike provides common data and functionality for OptFDNewton, OptQNewton, and OptNewton methods.

Note:
Modified by P.J. Williams, pwillia@sandia.gov
Date:
Last modified 03/2007


Constructor & Destructor Documentation

OPTPP::OptNewtonLike::OptNewtonLike (  )  [inline]

Default Constructor

See also:
OptNewtonLike(int n)

OptNewtonLike(int n, UPDATEFCN u)

OptNewtonLike(int n, TOLS t)

OPTPP::OptNewtonLike::OptNewtonLike ( int  n  )  [inline]

Parameters:
n an integer argument.

OPTPP::OptNewtonLike::OptNewtonLike ( int  n,
UPDATEFCN  u 
) [inline]

Parameters:
n an integer argument.
u a function pointer.

OPTPP::OptNewtonLike::OptNewtonLike ( int  n,
TOLS  t 
) [inline]

Parameters:
n an integer argument.
t tolerance class reference.

virtual OPTPP::OptNewtonLike::~OptNewtonLike (  )  [inline, virtual]

Destructor


Member Function Documentation

virtual void OPTPP::OptNewtonLike::acceptStep ( int  k,
int  step_type 
) [inline, virtual]

Note:
Pure virtual functions

Each derived class must define these functions for themselves

Implements OPTPP::OptimizeClass.

int OPTPP::OptNewtonLike::checkAnalyticFDGrad (  ) 

Compare the analytic gradient with the finite difference gradient.

int OPTPP::OptNewtonLike::checkConvg (  )  [virtual]

Check to see if algorithm satisfies the convergence criterion.

Implements OPTPP::OptimizeClass.

int OPTPP::OptNewtonLike::checkDeriv (  )  [virtual]

Compare the analytic gradient with the finite difference gradient.

Reimplemented in OPTPP::OptFDNewton, OPTPP::OptNewton, and OPTPP::OptQNewton.

virtual NEWMAT::ColumnVector OPTPP::OptNewtonLike::computeSearch ( NEWMAT::SymmetricMatrix &  H  )  [inline, virtual]

Compute Newton direction.

Implements OPTPP::OptimizeClass.

virtual int OPTPP::OptNewtonLike::computeStep ( NEWMAT::ColumnVector  sk  )  [virtual]

Compute the step length along the Newton direction.

void OPTPP::OptNewtonLike::defaultAcceptStep ( int  iter,
int  step_type 
) [protected, virtual]

Provide default implementation of AcceptStep

Reimplemented from OPTPP::OptimizeClass.

NEWMAT::ColumnVector OPTPP::OptNewtonLike::defaultComputeSearch ( NEWMAT::SymmetricMatrix &   )  [protected]

DerivOption OPTPP::OptNewtonLike::getDerivOption (  )  const [inline]

Set the type of finite difference routine.

int OPTPP::OptNewtonLike::getFevals (  )  const [inline]

Returns:
The number of function evaluations

int OPTPP::OptNewtonLike::getGevals (  )  const [inline]

Returns:
The number of gradient evaluations

real OPTPP::OptNewtonLike::getGradMult (  )  const [inline]

Returns:
Gradient multiplier to compute TR_size

NEWMAT::SymmetricMatrix OPTPP::OptNewtonLike::getHessian (  )  const [inline]

Returns:
Hessian matrix

int OPTPP::OptNewtonLike::getSearchSize (  )  const [inline]

Returns:
Number of points in search scheme which is used in trustpds search strategy

SearchStrategy OPTPP::OptNewtonLike::getSearchStrategy (  )  const [inline]

Set globalization strategy for optimization algorithms.

real OPTPP::OptNewtonLike::getTRSize (  )  const [inline]

Returns:
Radius of trust-region

bool OPTPP::OptNewtonLike::getWarmStart (  )  const [inline]

void OPTPP::OptNewtonLike::initHessian (  )  [virtual]

Compute the Hessian or its approximation at the initial point.

Reimplemented in OPTPP::OptNewton.

void OPTPP::OptNewtonLike::initOpt (  )  [virtual]

Initialize algorithmic parameters.

double OPTPP::OptNewtonLike::initTrustRegionSize (  )  const [virtual]

the trustregion or trustpds globalization strategies are selected

virtual NLP1* OPTPP::OptNewtonLike::nlprob (  )  const [protected, pure virtual]

Returns pointer to an NLP1 object.

Implemented in OPTPP::OptNewton1Deriv, and OPTPP::OptNewton2Deriv.

void OPTPP::OptNewtonLike::optimize (  )  [virtual]

Invoke Newton's method on an unconstrained problem.

Implements OPTPP::OptimizeClass.

void OPTPP::OptNewtonLike::printStatus ( char *   )  [virtual]

Print status of unconstrained Newton's method.

Implements OPTPP::OptimizeClass.

Reimplemented in OPTPP::OptNewton.

void OPTPP::OptNewtonLike::readOptInput (  )  [virtual]

Read user-specified input options from a file.

Implements OPTPP::OptimizeClass.

void OPTPP::OptNewtonLike::reset (  )  [virtual]

Reset parameters.

Implements OPTPP::OptimizeClass.

void OPTPP::OptNewtonLike::setDerivOption ( DerivOption  d  )  [inline]

Returns:
Type of finite difference approximation.

void OPTPP::OptNewtonLike::setGradMult ( real  tau  )  [inline]

Set gradient multiplier which is used to compute trust-region radius.

void OPTPP::OptNewtonLike::setHessian ( NEWMAT::SymmetricMatrix &  H  )  [inline]

Store the current Hessian matrix.

void OPTPP::OptNewtonLike::setSearchSize ( int  sss  )  [inline]

Set number of points in search scheme for trust-pds search strategy.

void OPTPP::OptNewtonLike::setSearchStrategy ( SearchStrategy  s  )  [inline]

Returns:
Globalization strategy for optimization algorithms

void OPTPP::OptNewtonLike::setTRSize ( real  delta  )  [inline]

Set trust-region radius.

virtual NEWMAT::SymmetricMatrix OPTPP::OptNewtonLike::updateH ( NEWMAT::SymmetricMatrix &  H,
int  k 
) [pure virtual]

Compute the Hessian of the objective function or its approximation at the current point.

Implemented in OPTPP::OptFDNewton, OPTPP::OptNewton, and OPTPP::OptQNewton.

virtual void OPTPP::OptNewtonLike::updateModel ( int  k,
int  ndim,
NEWMAT::ColumnVector  x 
) [inline, virtual]

Implements OPTPP::OptimizeClass.

void OPTPP::OptNewtonLike::UseWarmStart ( NEWMAT::SymmetricMatrix &  H  )  [inline]


Friends And Related Function Documentation

int trustpds ( NLP1 ,
ostream *  ,
NEWMAT::SymmetricMatrix &  ,
NEWMAT::ColumnVector &  ,
NEWMAT::ColumnVector &  ,
real ,
real ,
real  stpmax,
real  stpmin,
int   
) [friend]

int trustregion ( NLP1 ,
ostream *  ,
NEWMAT::SymmetricMatrix &  ,
NEWMAT::ColumnVector &  ,
NEWMAT::ColumnVector &  ,
real ,
real ,
real  stpmax,
real  stpmin 
) [friend]


Member Data Documentation

DerivOption OPTPP::OptNewtonLike::finitediff [protected]

User-specified derivative option.

NEWMAT::ColumnVector OPTPP::OptNewtonLike::gprev [protected]

Gradient at the prev. iteration.

int OPTPP::OptNewtonLike::grad_evals [protected]

Number of gradient evaluations.

real OPTPP::OptNewtonLike::gradMult [protected]

Gradient multiplier to compute TR_size.

NEWMAT::SymmetricMatrix OPTPP::OptNewtonLike::Hessian [protected]

Current Hessian.

int OPTPP::OptNewtonLike::searchSize [protected]

Search pattern size for TRPDS.

SearchStrategy OPTPP::OptNewtonLike::strategy [protected]

User-specified globalization Strategy.

real OPTPP::OptNewtonLike::TR_size [protected]

Trust region radius.

bool OPTPP::OptNewtonLike::WarmStart [protected]


The documentation for this class was generated from the following files:
Bug Reports    OPT++ Developers    Copyright Information    GNU Lesser General Public License
Documentation, generated by , last revised August 30, 2006.