ROL
Classes | Public Member Functions | Private 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  BuildC
 
class  isGreater
 
class  LowerBinding
 
class  PruneBinding
 
class  SetZeroEntry
 
class  UpperBinding
 

Public Member Functions

 Bounds (const Vector< Real > &x, bool isLower=true, Real scale=1, Real feasTol=std::sqrt(ROL_EPSILON< Real >()))
 
 Bounds (const Ptr< Vector< Real >> &x_lo, const Ptr< Vector< Real >> &x_up, const Real scale=1, const Real feasTol=std::sqrt(ROL_EPSILON< Real >()))
 
void project (Vector< Real > &x) override
 Project optimization variables onto the bounds. More...
 
void projectInterior (Vector< Real > &x) override
 Project optimization variables into the interior of the feasible set. More...
 
void pruneUpperActive (Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0)) override
 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 xeps=Real(0), Real geps=Real(0)) override
 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=Real(0)) override
 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 xeps=Real(0), Real geps=Real(0)) override
 Set variables to zero if they correspond to the \(\epsilon\)-binding set. More...
 
bool isFeasible (const Vector< Real > &v) override
 Check if the vector, v, is feasible. More...
 
void applyInverseScalingFunction (Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const override
 Apply inverse scaling function. More...
 
void applyScalingFunctionJacobian (Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const override
 Apply scaling function Jacobian. More...
 
- Public Member Functions inherited from ROL::BoundConstraint< Real >
virtual ~BoundConstraint ()
 
 BoundConstraint (void)
 
 BoundConstraint (const Vector< Real > &x)
 
virtual const Ptr< const
Vector< Real > > 
getLowerBound (void) const
 Return the ref count pointer to the lower bound vector. More...
 
virtual const Ptr< const
Vector< Real > > 
getUpperBound (void) const
 Return the ref count pointer to the upper bound vector. 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=Real(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 xeps=Real(0), Real geps=Real(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=Real(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=Real(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 xeps=Real(0), Real geps=Real(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 xeps=Real(0), Real geps=Real(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=Real(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 xeps=Real(0), Real geps=Real(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 Member Functions

void buildScalingFunction (Vector< Real > &d, const Vector< Real > &x, const Vector< Real > &g) const
 

Private Attributes

const Real scale_
 
const Real feasTol_
 
Ptr< Vector< Real > > mask_
 
Real min_diff_
 
Elementwise::ReductionMin< Real > minimum_
 
Elementwise::ReductionMax< Real > maximum_
 
ROL::Bounds::isGreater isGreater_
 
ROL::Bounds::PruneBinding prune_
 
ROL::Bounds::BuildC buildC_
 
ROL::Bounds::SetZeroEntry setZeroEntry_
 

Additional Inherited Members

- Protected Member Functions inherited from ROL::BoundConstraint< Real >
Real computeInf (const Vector< Real > &x) const
 
- Protected Attributes inherited from ROL::BoundConstraint< Real >
Ptr< Vector< Real > > lower_
 
Ptr< Vector< Real > > upper_
 

Detailed Description

template<typename 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<typename Real >
Bounds::Bounds ( const Vector< Real > &  x,
bool  isLower = true,
Real  scale = 1,
Real  feasTol = std::sqrt(ROL_EPSILON<Real>()) 
)
template<typename Real >
Bounds::Bounds ( const Ptr< Vector< Real >> &  x_lo,
const Ptr< Vector< Real >> &  x_up,
const Real  scale = 1,
const Real  feasTol = std::sqrt(ROL_EPSILON<Real>()) 
)

Member Function Documentation

template<typename Real >
void Bounds::buildScalingFunction ( Vector< Real > &  d,
const Vector< Real > &  x,
const Vector< Real > &  g 
) const
private
template<typename Real >
void Bounds::project ( Vector< Real > &  x)
overridevirtual

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 91 of file ROL_Bounds_Def.hpp.

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

template<typename Real >
void Bounds::projectInterior ( Vector< Real > &  x)
overridevirtual

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 109 of file ROL_Bounds_Def.hpp.

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

Referenced by testRandomInputs().

template<typename Real >
void Bounds::pruneUpperActive ( Vector< Real > &  v,
const Vector< Real > &  x,
Real  eps = Real(0) 
)
overridevirtual

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) \ge 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 156 of file ROL_Bounds_Def.hpp.

References ROL::Vector< Real >::applyBinary().

template<typename Real >
void Bounds::pruneUpperActive ( Vector< Real > &  v,
const Vector< Real > &  g,
const Vector< Real > &  x,
Real  xeps = Real(0),
Real  geps = Real(0) 
)
overridevirtual

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) \ge b(\xi)-\epsilon_x,\; g(\xi) < -\epsilon_g \,\}. \]

Parameters
[out]vis the variable to be pruned.
[in]gis the negative search direction.
[in]xis the current optimization variable.
[in]xepsis the active-set tolerance \(\epsilon_x\).
[in]gepsis the binding-set tolerance \(\epsilon_g\).

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 169 of file ROL_Bounds_Def.hpp.

References ROL::Vector< Real >::applyBinary().

template<typename Real >
void Bounds::pruneLowerActive ( Vector< Real > &  v,
const Vector< Real > &  x,
Real  eps = Real(0) 
)
overridevirtual

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) \le 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 184 of file ROL_Bounds_Def.hpp.

References ROL::Vector< Real >::applyBinary().

template<typename Real >
void Bounds::pruneLowerActive ( Vector< Real > &  v,
const Vector< Real > &  g,
const Vector< Real > &  x,
Real  xeps = Real(0),
Real  geps = Real(0) 
)
overridevirtual

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) \le 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]xepsis the active-set tolerance \(\epsilon_x\).
[in]gepsis the binding-set tolerance \(\epsilon_g\).

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 197 of file ROL_Bounds_Def.hpp.

References ROL::Vector< Real >::applyBinary().

template<typename Real >
bool Bounds::isFeasible ( const Vector< Real > &  v)
overridevirtual

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 212 of file ROL_Bounds_Def.hpp.

template<typename Real >
void Bounds::applyInverseScalingFunction ( Vector< Real > &  dv,
const Vector< Real > &  v,
const Vector< Real > &  x,
const Vector< Real > &  g 
) const
overridevirtual

Apply inverse scaling function.

This function applies the inverse scaling function \(d(x,g)\) to a vector \(v\), i.e., the output is \(\mathrm{diag}(d(x,g)^{-1})v\). The scaling function must satisfy: (i) \(d(x,g)_i = 0\) if \(x_i = a_i\) and \(g_i \ge 0\); (ii) \(d(x,g)_i = 0\) if \(x_i = b_i\) and \(g_i \le 0\); and (iii) \(d(x,g)_i > 0\) otherwise.

Parameters
[out]dvis the inverse scaling function applied to v.
[in]vis the vector being scaled.
[in]xis the primal vector at which the scaling function is evaluated.
[in]gis the dual vector at which the scaling function is evaluated.

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 289 of file ROL_Bounds_Def.hpp.

References ROL::Vector< Real >::applyBinary(), and ROL::Bounds< Real >::buildScalingFunction().

Referenced by testRandomInputs().

template<typename Real >
void Bounds::applyScalingFunctionJacobian ( Vector< Real > &  dv,
const Vector< Real > &  v,
const Vector< Real > &  x,
const Vector< Real > &  g 
) const
overridevirtual

Apply scaling function Jacobian.

This function applies the Jacobian of the scaling function \(d(x,g)\) to a vector \(v\). The output is \(\mathrm{diag}(d_x(x,g)g)v\). The scaling function must satisfy: (i) \(d(x,g)_i = 0\) if \(x_i = a_i\) and \(g_i \ge 0\); (ii) \(d(x,g)_i = 0\) if \(x_i = b_i\) and \(g_i \le 0\); and (iii) \(d(x,g)_i > 0\) otherwise.

Parameters
[out]dvis the scaling function Jacobian applied to v.
[in]vis the vector being scaled.
[in]xis the primal vector at which the scaling function is evaluated.
[in]gis the dual vector at which the scaling function is evaluated.

Reimplemented from ROL::BoundConstraint< Real >.

Definition at line 295 of file ROL_Bounds_Def.hpp.

References ROL::Vector< Real >::applyBinary(), ROL::Vector< Real >::axpy(), ROL::Bounds< Real >::buildScalingFunction(), ROL::Vector< Real >::plus(), and ROL::Vector< Real >::setScalar().

Referenced by testRandomInputs().

Member Data Documentation

template<typename Real>
const Real ROL::Bounds< Real >::scale_
private

Definition at line 61 of file ROL_Bounds.hpp.

template<typename Real>
const Real ROL::Bounds< Real >::feasTol_
private

Definition at line 62 of file ROL_Bounds.hpp.

template<typename Real>
Ptr<Vector<Real> > ROL::Bounds< Real >::mask_
private

Definition at line 67 of file ROL_Bounds.hpp.

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

template<typename Real>
Real ROL::Bounds< Real >::min_diff_
private

Definition at line 69 of file ROL_Bounds.hpp.

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

template<typename Real>
Elementwise::ReductionMin<Real> ROL::Bounds< Real >::minimum_
private

Definition at line 71 of file ROL_Bounds.hpp.

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

template<typename Real>
Elementwise::ReductionMax<Real> ROL::Bounds< Real >::maximum_
private

Definition at line 72 of file ROL_Bounds.hpp.

template<typename Real>
ROL::Bounds::isGreater ROL::Bounds< Real >::isGreater_
private
template<typename Real>
ROL::Bounds::PruneBinding ROL::Bounds< Real >::prune_
private
template<typename Real>
ROL::Bounds::BuildC ROL::Bounds< Real >::buildC_
private
template<typename Real>
ROL::Bounds::SetZeroEntry ROL::Bounds< Real >::setZeroEntry_
private

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