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

#include <ROL_RiskBoundConstraint.hpp>

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

Public Member Functions

 RiskBoundConstraint (Ptr< ParameterList > &parlist, const Ptr< BoundConstraint< Real >> &bc=nullPtr)
 
 RiskBoundConstraint (std::vector< Ptr< ParameterList >> &parlist, const Ptr< BoundConstraint< Real >> &bc=nullPtr)
 
 RiskBoundConstraint (Ptr< ParameterList > &parlistObj, std::vector< Ptr< ParameterList >> &parlistCon, const Ptr< BoundConstraint< Real >> &bc=nullPtr)
 
 RiskBoundConstraint (const Ptr< BoundConstraint< Real >> &bc)
 
void project (Vector< Real > &x)
 Project optimization variables onto the bounds. More...
 
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=Real(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 xeps=Real(0), Real geps=Real(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=Real(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 xeps=Real(0), Real geps=Real(0))
 Set variables to zero if they correspond to the \(\epsilon\)-binding set. More...
 
const Ptr< const Vector< Real > > getLowerBound (void) const
 Return the ref count pointer to the lower bound vector. More...
 
const Ptr< const Vector< Real > > getUpperBound (void) const
 Return the ref count pointer to the upper bound vector. More...
 
bool isFeasible (const Vector< Real > &v)
 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
 Apply inverse scaling function. More...
 
void applyScalingFunctionJacobian (Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const
 Apply scaling function Jacobian. More...
 
- Public Member Functions inherited from ROL::BoundConstraint< Real >
virtual ~BoundConstraint ()
 
 BoundConstraint (void)
 
 BoundConstraint (const Vector< Real > &x)
 
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 setBoundInfo (ParameterList &parlist, int &nStat, std::vector< Real > &lower, std::vector< Real > &upper, bool &augmented, bool &activated)
 
bool buildObjStatBnd (Ptr< ParameterList > &parlist)
 
bool buildConStatBnd (std::vector< Ptr< ParameterList >> &parlist)
 

Private Attributes

Ptr< BoundConstraint< Real > > bc_
 
Ptr< StdBoundConstraint< Real > > statObj_bc_
 
std::vector< Real > lowerObj_
 
std::vector< Real > upperObj_
 
std::vector< Ptr
< StdBoundConstraint< Real > > > 
statCon_bc_
 
std::vector< std::vector< Real > > lowerCon_
 
std::vector< std::vector< Real > > upperCon_
 
bool augmentedObj_
 
bool activatedObj_
 
int nStatObj_
 
bool augmentedCon_
 
std::vector< bool > activatedCon_
 
std::vector< int > nStatCon_
 
bool isLOinitialized_
 
bool isHIinitialized_
 
Ptr< RiskVector< Real > > lo_
 
Ptr< RiskVector< Real > > hi_
 

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<class Real>
class ROL::RiskBoundConstraint< Real >

Definition at line 54 of file ROL_RiskBoundConstraint.hpp.

Constructor & Destructor Documentation

template<class Real >
ROL::RiskBoundConstraint< Real >::RiskBoundConstraint ( Ptr< ParameterList > &  parlist,
const Ptr< BoundConstraint< Real >> &  bc = nullPtr 
)
inline
template<class Real >
ROL::RiskBoundConstraint< Real >::RiskBoundConstraint ( std::vector< Ptr< ParameterList >> &  parlist,
const Ptr< BoundConstraint< Real >> &  bc = nullPtr 
)
inline
template<class Real >
ROL::RiskBoundConstraint< Real >::RiskBoundConstraint ( Ptr< ParameterList > &  parlistObj,
std::vector< Ptr< ParameterList >> &  parlistCon,
const Ptr< BoundConstraint< Real >> &  bc = nullPtr 
)
inline
template<class Real >
ROL::RiskBoundConstraint< Real >::RiskBoundConstraint ( const Ptr< BoundConstraint< Real >> &  bc)
inline

Member Function Documentation

template<class Real >
void ROL::RiskBoundConstraint< Real >::setBoundInfo ( ParameterList &  parlist,
int &  nStat,
std::vector< Real > &  lower,
std::vector< Real > &  upper,
bool &  augmented,
bool &  activated 
)
inlineprivate
template<class Real >
bool ROL::RiskBoundConstraint< Real >::buildObjStatBnd ( Ptr< ParameterList > &  parlist)
inlineprivate
template<class Real >
bool ROL::RiskBoundConstraint< Real >::buildConStatBnd ( std::vector< Ptr< ParameterList >> &  parlist)
inlineprivate
template<class Real >
void ROL::RiskBoundConstraint< 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 236 of file ROL_RiskBoundConstraint.hpp.

References ROL::RiskBoundConstraint< Real >::activatedCon_, ROL::RiskBoundConstraint< Real >::activatedObj_, ROL::RiskBoundConstraint< Real >::augmentedCon_, ROL::RiskBoundConstraint< Real >::augmentedObj_, ROL::RiskBoundConstraint< Real >::bc_, ROL::RiskBoundConstraint< Real >::statCon_bc_, and ROL::RiskBoundConstraint< Real >::statObj_bc_.

template<class Real >
void ROL::RiskBoundConstraint< 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 256 of file ROL_RiskBoundConstraint.hpp.

References ROL::RiskBoundConstraint< Real >::activatedCon_, ROL::RiskBoundConstraint< Real >::activatedObj_, ROL::RiskBoundConstraint< Real >::augmentedCon_, ROL::RiskBoundConstraint< Real >::augmentedObj_, ROL::RiskBoundConstraint< Real >::bc_, ROL::RiskBoundConstraint< Real >::statCon_bc_, and ROL::RiskBoundConstraint< Real >::statObj_bc_.

template<class Real >
void ROL::RiskBoundConstraint< Real >::pruneUpperActive ( Vector< Real > &  v,
const Vector< Real > &  x,
Real  eps = Real(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) \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 276 of file ROL_RiskBoundConstraint.hpp.

References ROL::RiskBoundConstraint< Real >::activatedCon_, ROL::RiskBoundConstraint< Real >::activatedObj_, ROL::RiskBoundConstraint< Real >::augmentedCon_, ROL::RiskBoundConstraint< Real >::augmentedObj_, ROL::RiskBoundConstraint< Real >::bc_, ROL::RiskBoundConstraint< Real >::statCon_bc_, and ROL::RiskBoundConstraint< Real >::statObj_bc_.

template<class Real >
void ROL::RiskBoundConstraint< Real >::pruneUpperActive ( Vector< Real > &  v,
const Vector< Real > &  g,
const Vector< Real > &  x,
Real  xeps = Real(0),
Real  geps = Real(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) \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 299 of file ROL_RiskBoundConstraint.hpp.

References ROL::RiskBoundConstraint< Real >::activatedCon_, ROL::RiskBoundConstraint< Real >::activatedObj_, ROL::RiskBoundConstraint< Real >::augmentedCon_, ROL::RiskBoundConstraint< Real >::augmentedObj_, ROL::RiskBoundConstraint< Real >::bc_, ROL::RiskBoundConstraint< Real >::statCon_bc_, and ROL::RiskBoundConstraint< Real >::statObj_bc_.

template<class Real >
void ROL::RiskBoundConstraint< Real >::pruneLowerActive ( Vector< Real > &  v,
const Vector< Real > &  x,
Real  eps = Real(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) \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 325 of file ROL_RiskBoundConstraint.hpp.

References ROL::RiskBoundConstraint< Real >::activatedCon_, ROL::RiskBoundConstraint< Real >::activatedObj_, ROL::RiskBoundConstraint< Real >::augmentedCon_, ROL::RiskBoundConstraint< Real >::augmentedObj_, ROL::RiskBoundConstraint< Real >::bc_, ROL::RiskBoundConstraint< Real >::statCon_bc_, and ROL::RiskBoundConstraint< Real >::statObj_bc_.

template<class Real >
void ROL::RiskBoundConstraint< Real >::pruneLowerActive ( Vector< Real > &  v,
const Vector< Real > &  g,
const Vector< Real > &  x,
Real  xeps = Real(0),
Real  geps = Real(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) \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 348 of file ROL_RiskBoundConstraint.hpp.

References ROL::RiskBoundConstraint< Real >::activatedCon_, ROL::RiskBoundConstraint< Real >::activatedObj_, ROL::RiskBoundConstraint< Real >::augmentedCon_, ROL::RiskBoundConstraint< Real >::augmentedObj_, ROL::RiskBoundConstraint< Real >::bc_, ROL::RiskBoundConstraint< Real >::statCon_bc_, and ROL::RiskBoundConstraint< Real >::statObj_bc_.

template<class Real >
const Ptr<const Vector<Real> > ROL::RiskBoundConstraint< Real >::getLowerBound ( void  ) const
inlinevirtual
template<class Real >
const Ptr<const Vector<Real> > ROL::RiskBoundConstraint< Real >::getUpperBound ( void  ) const
inlinevirtual
template<class Real >
bool ROL::RiskBoundConstraint< Real >::isFeasible ( const Vector< Real > &  v)
inlinevirtual
template<class Real >
void ROL::RiskBoundConstraint< Real >::applyInverseScalingFunction ( Vector< Real > &  dv,
const Vector< Real > &  v,
const Vector< Real > &  x,
const Vector< Real > &  g 
) const
inlinevirtual

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 428 of file ROL_RiskBoundConstraint.hpp.

References ROL::RiskBoundConstraint< Real >::activatedCon_, ROL::RiskBoundConstraint< Real >::activatedObj_, ROL::RiskBoundConstraint< Real >::augmentedCon_, ROL::RiskBoundConstraint< Real >::augmentedObj_, ROL::RiskBoundConstraint< Real >::bc_, ROL::RiskBoundConstraint< Real >::statCon_bc_, and ROL::RiskBoundConstraint< Real >::statObj_bc_.

template<class Real >
void ROL::RiskBoundConstraint< Real >::applyScalingFunctionJacobian ( Vector< Real > &  dv,
const Vector< Real > &  v,
const Vector< Real > &  x,
const Vector< Real > &  g 
) const
inlinevirtual

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 457 of file ROL_RiskBoundConstraint.hpp.

References ROL::RiskBoundConstraint< Real >::activatedCon_, ROL::RiskBoundConstraint< Real >::activatedObj_, ROL::RiskBoundConstraint< Real >::augmentedCon_, ROL::RiskBoundConstraint< Real >::augmentedObj_, ROL::RiskBoundConstraint< Real >::bc_, ROL::RiskBoundConstraint< Real >::statCon_bc_, and ROL::RiskBoundConstraint< Real >::statObj_bc_.

Member Data Documentation

template<class Real >
Ptr<BoundConstraint<Real> > ROL::RiskBoundConstraint< Real >::bc_
private
template<class Real >
Ptr<StdBoundConstraint<Real> > ROL::RiskBoundConstraint< Real >::statObj_bc_
private
template<class Real >
std::vector<Real> ROL::RiskBoundConstraint< Real >::lowerObj_
private
template<class Real >
std::vector<Real> ROL::RiskBoundConstraint< Real >::upperObj_
private
template<class Real >
std::vector<Ptr<StdBoundConstraint<Real> > > ROL::RiskBoundConstraint< Real >::statCon_bc_
private
template<class Real >
std::vector<std::vector<Real> > ROL::RiskBoundConstraint< Real >::lowerCon_
private
template<class Real >
std::vector<std::vector<Real> > ROL::RiskBoundConstraint< Real >::upperCon_
private
template<class Real >
bool ROL::RiskBoundConstraint< Real >::augmentedObj_
private
template<class Real >
bool ROL::RiskBoundConstraint< Real >::activatedObj_
private
template<class Real >
int ROL::RiskBoundConstraint< Real >::nStatObj_
private
template<class Real >
bool ROL::RiskBoundConstraint< Real >::augmentedCon_
private
template<class Real >
std::vector<bool> ROL::RiskBoundConstraint< Real >::activatedCon_
private
template<class Real >
std::vector<int> ROL::RiskBoundConstraint< Real >::nStatCon_
private
template<class Real >
bool ROL::RiskBoundConstraint< Real >::isLOinitialized_
mutableprivate
template<class Real >
bool ROL::RiskBoundConstraint< Real >::isHIinitialized_
mutableprivate
template<class Real >
Ptr<RiskVector<Real> > ROL::RiskBoundConstraint< Real >::lo_
mutableprivate
template<class Real >
Ptr<RiskVector<Real> > ROL::RiskBoundConstraint< Real >::hi_
mutableprivate

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