ROL
ROL_StdBoundConstraint.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
15#ifndef ROL_STDBOUNDCONSTRAINT_HPP
16#define ROL_STDBOUNDCONSTRAINT_HPP
17
18#include "ROL_StdVector.hpp"
20
21namespace ROL {
22
23template<class Real>
24class StdBoundConstraint : public BoundConstraint<Real> {
25private:
26 int dim_;
27
28 std::vector<Real> x_lo_;
29 std::vector<Real> x_up_;
30
31 const Real scale_;
32 const Real feasTol_;
33
34 using BoundConstraint<Real>::lower_;
35 using BoundConstraint<Real>::upper_;
36
38
39 inline Real buildC(int i) const {
40 const Real zeta(0.5), kappa(1);
41 return std::min(zeta*(x_up_[i] - x_lo_[i]), kappa);
42 }
43
44 inline Real sgn(Real x) const {
45 const Real zero(0), one(1);
46 return x > zero ? one : (x < zero ? -one : zero);
47 }
48
49 void buildScalingFunction(Vector<Real> &d, const Vector<Real> &x, const Vector<Real> &g) const;
50
51public:
52 StdBoundConstraint(std::vector<Real> &x,
53 bool isLower = false,
54 Real scale = Real(1),
55 const Real feasTol = std::sqrt(ROL_EPSILON<Real>()));
56
57 StdBoundConstraint(std::vector<Real> &l,
58 std::vector<Real> &u,
59 Real scale = Real(1),
60 const Real feasTol = std::sqrt(ROL_EPSILON<Real>()));
61
62 void project( Vector<Real> &x ) override;
63
64 void projectInterior( Vector<Real> &x ) override;
65
66 void pruneUpperActive(Vector<Real> &v, const Vector<Real> &x, Real eps = Real(0)) override;
67
68 void pruneUpperActive(Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real xeps = Real(0), Real geps = Real(0)) override;
69
70 void pruneLowerActive(Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real xeps = Real(0), Real geps = Real(0)) override;
71
72 void pruneLowerActive(Vector<Real> &v, const Vector<Real> &x, Real eps = Real(0)) override;
73
74 bool isFeasible( const Vector<Real> &v ) override;
75
76 void applyInverseScalingFunction( Vector<Real> &dv, const Vector<Real> &v, const Vector<Real> &x, const Vector<Real> &g) const override;
77
78 void applyScalingFunctionJacobian(Vector<Real> &dv, const Vector<Real> &v, const Vector<Real> &x, const Vector<Real> &g) const override;
79};
80
81}// End ROL Namespace
82
84
85#endif
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Contains definitions for std::vector bound constraints.
Provides the interface to apply upper and lower bound constraints.
Ptr< Vector< Real > > upper_
Ptr< Vector< Real > > lower_
void buildScalingFunction(Vector< Real > &d, const Vector< Real > &x, const Vector< Real > &g) const
bool isFeasible(const Vector< Real > &v) override
Check if the vector, v, is feasible.
void applyScalingFunctionJacobian(Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const override
Apply scaling function Jacobian.
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 -binding set.
void projectInterior(Vector< Real > &x) override
Project optimization variables into the interior of the feasible set.
void applyInverseScalingFunction(Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const override
Apply inverse scaling function.
void project(Vector< Real > &x) override
Project optimization variables onto the bounds.
void pruneUpperActive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0)) override
Set variables to zero if they correspond to the upper -active set.
Defines the linear algebra or vector space interface.