GlobiPack Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | List of all members
GlobiPack::LineSearchBase< Scalar > Class Template Referenceabstract

Base class for 1D linearsearch algorithms. More...

#include <GlobiPack_LineSearchBase.hpp>

Inheritance diagram for GlobiPack::LineSearchBase< Scalar >:
Inheritance graph
[legend]

Public Member Functions

virtual bool requiresBaseDeriv () const =0
 Determines if the linesearch algorithm requires the base derivative at Dphi(0) or not. More...
 
virtual bool requiresDerivEvals () const =0
 Determines if the linesearch algorithm requires that Dphi(alpha) can be computed or not. More...
 
virtual bool doLineSearch (const MeritFunc1DBase< Scalar > &phi, const PointEval1D< Scalar > &point_k, const Ptr< PointEval1D< Scalar > > &point_kp1, const Ptr< int > &numIters) const =0
 Called to perform a linesearch. More...
 
- Public Member Functions inherited from Teuchos::Describable
DescribableStreamManipulatorState describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default)
 
std::ostream & operator<< (std::ostream &os, const DescribableStreamManipulatorState &d)
 
virtual std::string description () const
 
virtual void describe (FancyOStream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 
void describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 
virtual ~Describable ()
 
 LabeledObject ()
 
virtual ~LabeledObject ()
 
virtual void setObjectLabel (const std::string &objectLabel)
 
virtual std::string getObjectLabel () const
 
- Public Member Functions inherited from Teuchos::VerboseObject< LineSearchBase< Scalar > >
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT
RCP< const ParameterList
getValidVerboseObjectSublist ()
 
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT
void 
setupVerboseObjectSublist (ParameterList *paramList)
 
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT
void 
readVerboseObjectSublist (ParameterList *paramList, RCP< FancyOStream > *oStream, EVerbosityLevel *verbLevel)
 
void readVerboseObjectSublist (ParameterList *paramList, VerboseObject< LineSearchBase< Scalar > > *verboseObject)
 
 VerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null)
 
virtual const VerboseObjectsetVerbLevel (const EVerbosityLevel verbLevel) const
 
virtual const VerboseObjectsetOverridingVerbLevel (const EVerbosityLevel verbLevel) const
 
virtual EVerbosityLevel getVerbLevel () const
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptor
virtual void setParameterList (const RCP< ParameterList > &paramList)=0
 
virtual RCP< ParameterListgetNonconstParameterList ()=0
 
virtual RCP< ParameterListunsetParameterList ()=0
 
virtual RCP< const ParameterListgetParameterList () const
 
virtual RCP< const ParameterListgetValidParameters () const
 

Additional Inherited Members

- Static Public Member Functions inherited from Teuchos::VerboseObject< LineSearchBase< Scalar > >
static void setDefaultVerbLevel (const EVerbosityLevel defaultVerbLevel)
 
static EVerbosityLevel getDefaultVerbLevel ()
 
- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default
 
- Protected Member Functions inherited from Teuchos::VerboseObject< LineSearchBase< Scalar > >
void initializeVerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null)
 

Detailed Description

template<typename Scalar>
class GlobiPack::LineSearchBase< Scalar >

Base class for 1D linearsearch algorithms.

ToDo: Finish Documentation!

Definition at line 62 of file GlobiPack_LineSearchBase.hpp.

Member Function Documentation

template<typename Scalar >
virtual bool GlobiPack::LineSearchBase< Scalar >::requiresBaseDeriv ( ) const
pure virtual

Determines if the linesearch algorithm requires the base derivative at Dphi(0) or not.

Implemented in GlobiPack::BrentsLineSearch< Scalar >, and GlobiPack::ArmijoPolyInterpLineSearch< Scalar >.

template<typename Scalar >
virtual bool GlobiPack::LineSearchBase< Scalar >::requiresDerivEvals ( ) const
pure virtual

Determines if the linesearch algorithm requires that Dphi(alpha) can be computed or not.

Implemented in GlobiPack::BrentsLineSearch< Scalar >, and GlobiPack::ArmijoPolyInterpLineSearch< Scalar >.

template<typename Scalar >
virtual bool GlobiPack::LineSearchBase< Scalar >::doLineSearch ( const MeritFunc1DBase< Scalar > &  phi,
const PointEval1D< Scalar > &  point_k,
const Ptr< PointEval1D< Scalar > > &  point_kp1,
const Ptr< int > &  numIters 
) const
pure virtual

Called to perform a linesearch.

Parameters
phi[in] The merit function object that will compute the merit function value phi(alpha) and/or derivative Dphi(alpha) at different points alpha. The last call to phi.eval(...) will always be at the value of point_kp1->alpha returned.
point_k[in] The evaluation of the merit function and optionally its derivative at alpha=0.0.
point_kp1[in/out] On input, point_kp1->alpha is the initial value to try out (usually 1.0 for most Newton-based algorithms). Also, point_kp1->phi must be computed at this value for alpha as well as point_kp1->Dphi if required. On output, point_kp1->alpha is the accepted value for a successful line search, or it will be the alpha for the minimum phi(alpha) found during a failed line search algorithm.

Preconditions:

point_k.alpha == 0.0

point_k.phi != PointEval1D<Scalar>::valNotGiven()

[this->requiresBaseDeriv()==true] point_k.Dphi != PointEval1D<Scalar>::valNotGiven()

[this->requiresBaseDeriv()==true] point_k.Dphi < 0.0 (throw Exceptions::NotDescentDirection)

[this->requiresDerivEvals()==true] phi.supportsDerivEvals()==true

!is_null(point_kp1)

point_kp1->phi != PointEval1D<Scalar>::valNotGiven()

[this->requiresDerivEvals()==true] point_kp1->Dphi != PointEval1D<Scalar>::valNotGiven()

Returns
true for successful line search or false for a line search failure.

This function computes the approximate minimum to 1D merit function phi(alpha). More specifically the following problem is approximately solved:

  min  phi(alpha)  s.t. alpha = [0, alpha_upper]<br>

For many lineserach algorithms, if the initial point_kp1->alpha satisfies the internally defined descent requirement, then it will typically be choosen over smaller values of point_kp1->alpha that may result in a greater reduction in the given merit function. Other linesearch implementions will actually seek an approximate minimizer.

If the maximum number of iterations is exceeded without finding an acceptable point, then the subclass object will return false and will return values of point_kp1->alpha and point_kp1->phi will be for the lowest value of phi_kp1 = phi(alpha_k) found. In this case, the last call to phi(alpha_k) will be this best value of phi_kp1.

Implemented in GlobiPack::BrentsLineSearch< Scalar >, and GlobiPack::ArmijoPolyInterpLineSearch< Scalar >.


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