ROL
Classes | Public Member Functions | Private Types | Private Attributes | List of all members
ROL::InteriorPointPenalty< Real > Class Template Reference

Provides the interface to evaluate the Interior Pointy log barrier penalty function with upper and lower bounds on some elements. More...

#include <ROL_InteriorPointPenalty.hpp>

+ Inheritance diagram for ROL::InteriorPointPenalty< Real >:

Classes

class  Mask
 
class  ModifiedDivide
 
class  ModifiedLogarithm
 
class  ModifiedReciprocal
 

Public Member Functions

 ~InteriorPointPenalty ()
 
 InteriorPointPenalty (const ROL::Ptr< Objective< Real > > &obj, const ROL::Ptr< BoundConstraint< Real > > &con, ROL::ParameterList &parlist)
 
Real getObjectiveValue (void) const
 
ROL::Ptr< Vector< Real > > getGradient (void) const
 
int getNumberFunctionEvaluations (void) const
 
int getNumberGradientEvaluations (void) const
 
ROL::Ptr< const Vector< Real > > getLowerMask (void) const
 
ROL::Ptr< const Vector< Real > > getUpperMask (void) const
 
void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update the interior point penalized objective. More...
 
Real value (const Vector< Real > &x, Real &tol)
 Compute value. More...
 
void gradient (Vector< Real > &g, const Vector< Real > &x, Real &tol)
 Compute gradient. More...
 
void hessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Compute action of Hessian on vector. More...
 
const ROL::Ptr< OBJgetObjective (void)
 
const ROL::Ptr< BNDgetBoundConstraint (void)
 
- Public Member Functions inherited from ROL::Objective< Real >
virtual ~Objective ()
 
virtual Real dirDeriv (const Vector< Real > &x, const Vector< Real > &d, Real &tol)
 Compute directional derivative. More...
 
virtual void invHessVec (Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply inverse Hessian approximation to vector. More...
 
virtual void precond (Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
 Apply preconditioner to vector. More...
 
virtual std::vector
< std::vector< Real > > 
checkGradient (const Vector< Real > &x, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference gradient check. More...
 
virtual std::vector
< std::vector< Real > > 
checkGradient (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &d, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference gradient check. More...
 
virtual std::vector
< std::vector< Real > > 
checkGradient (const Vector< Real > &x, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference gradient check with specified step sizes. More...
 
virtual std::vector
< std::vector< Real > > 
checkGradient (const Vector< Real > &x, const Vector< Real > &g, const Vector< Real > &d, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference gradient check with specified step sizes. More...
 
virtual std::vector
< std::vector< Real > > 
checkHessVec (const Vector< Real > &x, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference Hessian-applied-to-vector check. More...
 
virtual std::vector
< std::vector< Real > > 
checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const bool printToStream=true, std::ostream &outStream=std::cout, const int numSteps=ROL_NUM_CHECKDERIV_STEPS, const int order=1)
 Finite-difference Hessian-applied-to-vector check. More...
 
virtual std::vector
< std::vector< Real > > 
checkHessVec (const Vector< Real > &x, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference Hessian-applied-to-vector check with specified step sizes. More...
 
virtual std::vector
< std::vector< Real > > 
checkHessVec (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const std::vector< Real > &steps, const bool printToStream=true, std::ostream &outStream=std::cout, const int order=1)
 Finite-difference Hessian-applied-to-vector check with specified step sizes. More...
 
virtual std::vector< Real > checkHessSym (const Vector< Real > &x, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout)
 Hessian symmetry check. More...
 
virtual std::vector< Real > checkHessSym (const Vector< Real > &x, const Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &w, const bool printToStream=true, std::ostream &outStream=std::cout)
 Hessian symmetry check. More...
 
virtual void setParameter (const std::vector< Real > &param)
 

Private Types

typedef Vector< Real > V
 
typedef Objective< Real > OBJ
 
typedef BoundConstraint< Real > BND
 
typedef Elementwise::ValueSet
< Real > 
ValueSet
 

Private Attributes

const ROL::Ptr< OBJobj_
 
const ROL::Ptr< BNDbnd_
 
const ROL::Ptr< const Vlo_
 
const ROL::Ptr< const Vup_
 
ROL::Ptr< Vg_
 
ROL::Ptr< VmaskL_
 
ROL::Ptr< VmaskU_
 
ROL::Ptr< Va_
 
ROL::Ptr< Vb_
 
ROL::Ptr< Vc_
 
bool useLinearDamping_
 
Real mu_
 
Real kappaD_
 
Real fval_
 
int nfval_
 
int ngval_
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::Objective< Real >
const std::vector< Real > getParameter (void) const
 

Detailed Description

template<class Real>
class ROL::InteriorPointPenalty< Real >

Provides the interface to evaluate the Interior Pointy log barrier penalty function with upper and lower bounds on some elements.


Definition at line 63 of file ROL_InteriorPointPenalty.hpp.

Member Typedef Documentation

template<class Real >
typedef Vector<Real> ROL::InteriorPointPenalty< Real >::V
private

Definition at line 65 of file ROL_InteriorPointPenalty.hpp.

template<class Real >
typedef Objective<Real> ROL::InteriorPointPenalty< Real >::OBJ
private

Definition at line 66 of file ROL_InteriorPointPenalty.hpp.

template<class Real >
typedef BoundConstraint<Real> ROL::InteriorPointPenalty< Real >::BND
private

Definition at line 67 of file ROL_InteriorPointPenalty.hpp.

template<class Real >
typedef Elementwise::ValueSet<Real> ROL::InteriorPointPenalty< Real >::ValueSet
private

Definition at line 69 of file ROL_InteriorPointPenalty.hpp.

Constructor & Destructor Documentation

template<class Real >
ROL::InteriorPointPenalty< Real >::~InteriorPointPenalty ( )
inline

Definition at line 145 of file ROL_InteriorPointPenalty.hpp.

template<class Real >
ROL::InteriorPointPenalty< Real >::InteriorPointPenalty ( const ROL::Ptr< Objective< Real > > &  obj,
const ROL::Ptr< BoundConstraint< Real > > &  con,
ROL::ParameterList &  parlist 
)
inline

Member Function Documentation

template<class Real >
Real ROL::InteriorPointPenalty< Real >::getObjectiveValue ( void  ) const
inline
template<class Real >
ROL::Ptr<Vector<Real> > ROL::InteriorPointPenalty< Real >::getGradient ( void  ) const
inline
template<class Real >
int ROL::InteriorPointPenalty< Real >::getNumberFunctionEvaluations ( void  ) const
inline
template<class Real >
int ROL::InteriorPointPenalty< Real >::getNumberGradientEvaluations ( void  ) const
inline
template<class Real >
ROL::Ptr<const Vector<Real> > ROL::InteriorPointPenalty< Real >::getLowerMask ( void  ) const
inline
template<class Real >
ROL::Ptr<const Vector<Real> > ROL::InteriorPointPenalty< Real >::getUpperMask ( void  ) const
inline
template<class Real >
void ROL::InteriorPointPenalty< Real >::update ( const Vector< Real > &  x,
bool  flag = true,
int  iter = -1 
)
inlinevirtual

Update the interior point penalized objective.

This function updates the log barrier penalized function at new iterations.

Parameters
[in]xis the new iterate.
[in]flagis true if the iterate has changed.
[in]iteris the outer algorithm iterations count.

Reimplemented from ROL::Objective< Real >.

Definition at line 213 of file ROL_InteriorPointPenalty.hpp.

References ROL::InteriorPointPenalty< Real >::bnd_, and ROL::InteriorPointPenalty< Real >::obj_.

template<class Real >
Real ROL::InteriorPointPenalty< Real >::value ( const Vector< Real > &  x,
Real &  tol 
)
inlinevirtual

Compute value.

This function returns the log barrier penalized objective value.

\[ \varphi_\mu(x) = f(x) - \mu \sum\limits_{i\int I_L} \ln(x_i-l_i) - \mu \sum\limits_{i\in I_Y} \ln(u_i-x_i) \]

Where \( I_L=\{i:l_i>-\infty\} \) and \( I_U = \{i:u_i<\infty\}\)

Parameters
[in]xis the current iterate.
[in]tolis a tolerance for interior point penalty computation.

Implements ROL::Objective< Real >.

Definition at line 229 of file ROL_InteriorPointPenalty.hpp.

References ROL::InteriorPointPenalty< Real >::a_, ROL::InteriorPointPenalty< Real >::b_, ROL::InteriorPointPenalty< Real >::c_, ROL::InteriorPointPenalty< Real >::fval_, ROL::InteriorPointPenalty< Real >::kappaD_, ROL::InteriorPointPenalty< Real >::lo_, ROL::InteriorPointPenalty< Real >::maskL_, ROL::InteriorPointPenalty< Real >::maskU_, ROL::InteriorPointPenalty< Real >::mu_, ROL::InteriorPointPenalty< Real >::nfval_, ROL::InteriorPointPenalty< Real >::obj_, ROL::InteriorPointPenalty< Real >::up_, and ROL::InteriorPointPenalty< Real >::useLinearDamping_.

template<class Real >
void ROL::InteriorPointPenalty< Real >::gradient ( Vector< Real > &  g,
const Vector< Real > &  x,
Real &  tol 
)
inlinevirtual
template<class Real >
void ROL::InteriorPointPenalty< Real >::hessVec ( Vector< Real > &  hv,
const Vector< Real > &  v,
const Vector< Real > &  x,
Real &  tol 
)
inlinevirtual

Compute action of Hessian on vector.

This function returns the log barrier penalized Hessian acting on a given vector.

Parameters
[out]hvis the Hessian-vector product.
[in]vis the given vector.
[in]xis the current iterate.
[in]tolis a tolerance for inexact log barrier penalty computation.

Reimplemented from ROL::Objective< Real >.

Definition at line 344 of file ROL_InteriorPointPenalty.hpp.

References ROL::InteriorPointPenalty< Real >::a_, ROL::Vector< Real >::axpy(), ROL::InteriorPointPenalty< Real >::b_, ROL::InteriorPointPenalty< Real >::lo_, ROL::InteriorPointPenalty< Real >::maskL_, ROL::InteriorPointPenalty< Real >::maskU_, ROL::InteriorPointPenalty< Real >::mu_, ROL::InteriorPointPenalty< Real >::obj_, and ROL::InteriorPointPenalty< Real >::up_.

template<class Real >
const ROL::Ptr<OBJ> ROL::InteriorPointPenalty< Real >::getObjective ( void  )
inline
template<class Real >
const ROL::Ptr<BND> ROL::InteriorPointPenalty< Real >::getBoundConstraint ( void  )
inline

Member Data Documentation

template<class Real >
const ROL::Ptr<OBJ> ROL::InteriorPointPenalty< Real >::obj_
private
template<class Real >
const ROL::Ptr<BND> ROL::InteriorPointPenalty< Real >::bnd_
private
template<class Real >
const ROL::Ptr<const V> ROL::InteriorPointPenalty< Real >::lo_
private
template<class Real >
const ROL::Ptr<const V> ROL::InteriorPointPenalty< Real >::up_
private
template<class Real >
ROL::Ptr<V> ROL::InteriorPointPenalty< Real >::g_
private
template<class Real >
ROL::Ptr<V> ROL::InteriorPointPenalty< Real >::maskL_
private
template<class Real >
ROL::Ptr<V> ROL::InteriorPointPenalty< Real >::maskU_
private
template<class Real >
ROL::Ptr<V> ROL::InteriorPointPenalty< Real >::a_
private
template<class Real >
ROL::Ptr<V> ROL::InteriorPointPenalty< Real >::b_
private
template<class Real >
ROL::Ptr<V> ROL::InteriorPointPenalty< Real >::c_
private
template<class Real >
bool ROL::InteriorPointPenalty< Real >::useLinearDamping_
private
template<class Real >
Real ROL::InteriorPointPenalty< Real >::mu_
private
template<class Real >
Real ROL::InteriorPointPenalty< Real >::kappaD_
private
template<class Real >
Real ROL::InteriorPointPenalty< Real >::fval_
private
template<class Real >
int ROL::InteriorPointPenalty< Real >::nfval_
private
template<class Real >
int ROL::InteriorPointPenalty< Real >::ngval_
private

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