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

Provides the elementwise interface to apply upper and lower bound constraints. More...

#include <ROL_Bounds.hpp>

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

Classes

class  Active
 
class  LowerBinding
 
class  PruneBinding
 
class  UpperBinding
 

Public Member Functions

 Bounds (const Vector< Real > &x, bool isLower=true, Real scale=1, Real feasTol=1e-2)
 
 Bounds (const ROL::Ptr< Vector< Real >> &x_lo, const ROL::Ptr< Vector< Real >> &x_up, const Real scale=1, const Real feasTol=1e-2)
 
virtual void project (Vector< Real > &x)
 Project optimization variables onto the bounds. More...
 
virtual void projectInterior (Vector< Real > &x)
 Project optimization variables into the interior of the feasible set. More...
 
void pruneUpperActive (Vector< Real > &v, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the upper \(\epsilon\)-active set. More...
 
void pruneUpperActive (Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the upper \(\epsilon\)-binding set. More...
 
void pruneLowerActive (Vector< Real > &v, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the lower \(\epsilon\)-active set. More...
 
void pruneLowerActive (Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the \(\epsilon\)-binding set. More...
 
const ROL::Ptr< const Vector
< Real > > 
getLowerBound (void) const
 Return the ref count pointer to the lower bound vector. More...
 
const ROL::Ptr< const Vector
< Real > > 
getUpperBound (void) const
 Return the ref count pointer to the upper bound vector. More...
 
virtual bool isFeasible (const Vector< Real > &v)
 Check if the vector, v, is feasible. More...
 
- Public Member Functions inherited from ROL::BoundConstraint< Real >
virtual ~BoundConstraint ()
 
 BoundConstraint (void)
 
 BoundConstraint (const Vector< Real > &x)
 
virtual void update (const Vector< Real > &x, bool flag=true, int iter=-1)
 Update bounds. More...
 
void activateLower (void)
 Turn on lower bound. More...
 
void activateUpper (void)
 Turn on upper bound. More...
 
void activate (void)
 Turn on bounds. More...
 
void deactivateLower (void)
 Turn off lower bound. More...
 
void deactivateUpper (void)
 Turn off upper bound. More...
 
void deactivate (void)
 Turn off bounds. More...
 
bool isLowerActivated (void) const
 Check if lower bound are on. More...
 
bool isUpperActivated (void) const
 Check if upper bound are on. More...
 
bool isActivated (void) const
 Check if bounds are on. More...
 
void pruneActive (Vector< Real > &v, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the \(\epsilon\)-active set. More...
 
void pruneActive (Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the \(\epsilon\)-binding set. More...
 
void pruneLowerInactive (Vector< Real > &v, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the \(\epsilon\)-inactive set. More...
 
void pruneUpperInactive (Vector< Real > &v, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the \(\epsilon\)-inactive set. More...
 
void pruneLowerInactive (Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the \(\epsilon\)-nonbinding set. More...
 
void pruneUpperInactive (Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the \(\epsilon\)-nonbinding set. More...
 
void pruneInactive (Vector< Real > &v, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the \(\epsilon\)-inactive set. More...
 
void pruneInactive (Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x, Real eps=0)
 Set variables to zero if they correspond to the \(\epsilon\)-nonbinding set. More...
 
void computeProjectedGradient (Vector< Real > &g, const Vector< Real > &x)
 Compute projected gradient. More...
 
void computeProjectedStep (Vector< Real > &v, const Vector< Real > &x)
 Compute projected step. More...
 

Private Attributes

const ROL::Ptr< Vector< Real > > x_lo_
 
const ROL::Ptr< Vector< Real > > x_up_
 
const Real scale_
 
const Real feasTol_
 
ROL::Ptr< Vector< Real > > mask_
 
Real min_diff_
 
Elementwise::ReductionMin< Real > minimum_
 
ROL::Bounds::PruneBinding prune_
 

Detailed Description

template<class Real>
class ROL::Bounds< Real >

Provides the elementwise interface to apply upper and lower bound constraints.

Definition at line 59 of file ROL_Bounds.hpp.

Constructor & Destructor Documentation

template<class Real>
ROL::Bounds< Real >::Bounds ( const Vector< Real > &  x,
bool  isLower = true,
Real  scale = 1,
Real  feasTol = 1e-2 
)
inline
template<class Real>
ROL::Bounds< Real >::Bounds ( const ROL::Ptr< Vector< Real >> &  x_lo,
const ROL::Ptr< Vector< Real >> &  x_up,
const Real  scale = 1,
const Real  feasTol = 1e-2 
)
inline

Member Function Documentation

template<class Real>
virtual void ROL::Bounds< Real >::project ( Vector< Real > &  x)
inlinevirtual

Project optimization variables onto the bounds.

This function implements the projection of \(x\) onto the bounds, i.e.,

\[ (P_{[a,b]}(x))(\xi) = \min\{b(\xi),\max\{a(\xi),x(\xi)\}\} \quad \text{for almost every }\xi\in\Xi. \]

Parameters
[in,out]xis the optimization variable.

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 146 of file ROL_Bounds.hpp.

References ROL::apply, ROL::Vector< Real >::applyBinary(), ROL::Bounds< Real >::x_lo_, and ROL::Bounds< Real >::x_up_.

template<class Real>
virtual void ROL::Bounds< Real >::projectInterior ( Vector< Real > &  x)
inlinevirtual

Project optimization variables into the interior of the feasible set.

This function implements the projection of \(x\) into the interior of the feasible set, i.e.,

\[ (\bar{P}_{[a,b]}(x))(\xi) \in (a(\xi),b(\xi)) \quad \text{for almost every }\xi\in\Xi. \]

Parameters
[in,out]xis the optimization variable.

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 163 of file ROL_Bounds.hpp.

References ROL::apply, ROL::Vector< Real >::applyBinary(), ROL::Bounds< Real >::feasTol_, ROL::Bounds< Real >::min_diff_, ROL::Bounds< Real >::x_lo_, and ROL::Bounds< Real >::x_up_.

template<class Real>
void ROL::Bounds< Real >::pruneUpperActive ( Vector< Real > &  v,
const Vector< Real > &  x,
Real  eps = 0 
)
inlinevirtual

Set variables to zero if they correspond to the upper \(\epsilon\)-active set.

This function sets \(v(\xi)=0\) if \(\xi\in\mathcal{A}^+_\epsilon(x)\). Here, the upper \(\epsilon\)-active set is defined as

\[ \mathcal{A}^+_\epsilon(x) = \{\,\xi\in\Xi\,:\,x(\xi) = b(\xi)-\epsilon\,\}. \]

Parameters
[out]vis the variable to be pruned.
[in]xis the current optimization variable.
[in]epsis the active-set tolerance \(\epsilon\).

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 209 of file ROL_Bounds.hpp.

References ROL::Vector< Real >::applyBinary(), ROL::Bounds< Real >::mask_, ROL::Bounds< Real >::min_diff_, ROL::Bounds< Real >::scale_, and ROL::Bounds< Real >::x_up_.

template<class Real>
void ROL::Bounds< Real >::pruneUpperActive ( Vector< Real > &  v,
const Vector< Real > &  g,
const Vector< Real > &  x,
Real  eps = 0 
)
inlinevirtual

Set variables to zero if they correspond to the upper \(\epsilon\)-binding set.

This function sets \(v(\xi)=0\) if \(\xi\in\mathcal{B}^+_\epsilon(x)\). Here, the upper \(\epsilon\)-binding set is defined as

\[ \mathcal{B}^+_\epsilon(x) = \{\,\xi\in\Xi\,:\,x(\xi) = b(\xi)-\epsilon,\; g(\xi) < 0 \,\}. \]

Parameters
[out]vis the variable to be pruned.
[in]gis the negative search direction.
[in]xis the current optimization variable.
[in]epsis the active-set tolerance \(\epsilon\).

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 221 of file ROL_Bounds.hpp.

References ROL::Vector< Real >::applyBinary(), ROL::Bounds< Real >::mask_, ROL::Bounds< Real >::min_diff_, ROL::Bounds< Real >::prune_, ROL::Bounds< Real >::scale_, and ROL::Bounds< Real >::x_up_.

template<class Real>
void ROL::Bounds< Real >::pruneLowerActive ( Vector< Real > &  v,
const Vector< Real > &  x,
Real  eps = 0 
)
inlinevirtual

Set variables to zero if they correspond to the lower \(\epsilon\)-active set.

This function sets \(v(\xi)=0\) if \(\xi\in\mathcal{A}^-_\epsilon(x)\). Here, the lower \(\epsilon\)-active set is defined as

\[ \mathcal{A}^-_\epsilon(x) = \{\,\xi\in\Xi\,:\,x(\xi) = a(\xi)+\epsilon\,\}. \]

Parameters
[out]vis the variable to be pruned.
[in]xis the current optimization variable.
[in]epsis the active-set tolerance \(\epsilon\).

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 235 of file ROL_Bounds.hpp.

References ROL::Vector< Real >::applyBinary(), ROL::Bounds< Real >::mask_, ROL::Bounds< Real >::min_diff_, ROL::Bounds< Real >::scale_, and ROL::Bounds< Real >::x_lo_.

template<class Real>
void ROL::Bounds< Real >::pruneLowerActive ( Vector< Real > &  v,
const Vector< Real > &  g,
const Vector< Real > &  x,
Real  eps = 0 
)
inlinevirtual

Set variables to zero if they correspond to the \(\epsilon\)-binding set.

This function sets \(v(\xi)=0\) if \(\xi\in\mathcal{B}^-_\epsilon(x)\). Here, the lower \(\epsilon\)-binding set is defined as

\[ \mathcal{B}^-_\epsilon(x) = \{\,\xi\in\Xi\,:\,x(\xi) = a(\xi)+\epsilon,\; g(\xi) > 0 \,\}. \]

Parameters
[out]vis the variable to be pruned.
[in]gis the negative search direction.
[in]xis the current optimization variable.
[in]epsis the active-set tolerance \(\epsilon\).

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 247 of file ROL_Bounds.hpp.

References ROL::Vector< Real >::applyBinary(), ROL::Bounds< Real >::mask_, ROL::Bounds< Real >::min_diff_, ROL::Bounds< Real >::prune_, ROL::Bounds< Real >::scale_, and ROL::Bounds< Real >::x_lo_.

template<class Real>
const ROL::Ptr<const Vector<Real> > ROL::Bounds< Real >::getLowerBound ( void  ) const
inlinevirtual

Return the ref count pointer to the lower bound vector.

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 261 of file ROL_Bounds.hpp.

References ROL::Bounds< Real >::x_lo_.

template<class Real>
const ROL::Ptr<const Vector<Real> > ROL::Bounds< Real >::getUpperBound ( void  ) const
inlinevirtual

Return the ref count pointer to the upper bound vector.

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 265 of file ROL_Bounds.hpp.

References ROL::Bounds< Real >::x_up_.

template<class Real>
virtual bool ROL::Bounds< Real >::isFeasible ( const Vector< Real > &  v)
inlinevirtual

Check if the vector, v, is feasible.

This function returns true if \(v = P_{[a,b]}(v)\).

Parameters
[in]vis the vector to be checked.

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 269 of file ROL_Bounds.hpp.

References ROL::Bounds< Real >::mask_, ROL::Bounds< Real >::minimum_, ROL::Bounds< Real >::x_lo_, and ROL::Bounds< Real >::x_up_.

Member Data Documentation

template<class Real>
const ROL::Ptr<Vector<Real> > ROL::Bounds< Real >::x_lo_
private
template<class Real>
const ROL::Ptr<Vector<Real> > ROL::Bounds< Real >::x_up_
private
template<class Real>
const Real ROL::Bounds< Real >::scale_
private
template<class Real>
const Real ROL::Bounds< Real >::feasTol_
private

Definition at line 64 of file ROL_Bounds.hpp.

Referenced by ROL::Bounds< Real >::projectInterior().

template<class Real>
ROL::Ptr<Vector<Real> > ROL::Bounds< Real >::mask_
private
template<class Real>
Real ROL::Bounds< Real >::min_diff_
private
template<class Real>
Elementwise::ReductionMin<Real> ROL::Bounds< Real >::minimum_
private

Definition at line 70 of file ROL_Bounds.hpp.

Referenced by ROL::Bounds< Real >::Bounds(), and ROL::Bounds< Real >::isFeasible().

template<class Real>
ROL::Bounds::PruneBinding ROL::Bounds< Real >::prune_
private

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