ROL
ROL_BoundConstraint_SimOpt.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
10#ifndef ROL_BOUND_CONSTRAINT_SIMOPT_H
11#define ROL_BOUND_CONSTRAINT_SIMOPT_H
12
14#include "ROL_Vector_SimOpt.hpp"
15#include "ROL_Types.hpp"
16#include <iostream>
17
34
35
36namespace ROL {
37
38template <class Real>
40private:
41 const Ptr<BoundConstraint<Real>> bnd1_, bnd2_;
42
43public:
45
51 const Ptr<BoundConstraint<Real>> &bnd2)
52 : bnd1_(bnd1), bnd2_(bnd2) {
53 if ( bnd1_->isActivated() || bnd2_->isActivated() )
55 else
57 }
58
67 bool optBnd = true)
68 : bnd1_(optBnd ? makePtr<BoundConstraint<Real>>() : bnd),
69 bnd2_(optBnd ? bnd : makePtr<BoundConstraint<Real>>()) {
70 if ( bnd1_->isActivated() || bnd2_->isActivated() )
72 else
74 }
75
85 const Vector<Real> &x,
86 bool optBnd = true)
87 : bnd1_(optBnd ? makePtr<BoundConstraint<Real>>(x) : bnd),
88 bnd2_(optBnd ? bnd : makePtr<BoundConstraint<Real>>(x)) {
89 if ( bnd1_->isActivated() || bnd2_->isActivated() )
91 else
93 }
94
104 Vector_SimOpt<Real> &xs = dynamic_cast<Vector_SimOpt<Real>&>(x);
105 if (bnd1_->isActivated()) bnd1_->project(*(xs.get_1()));
106 if (bnd2_->isActivated()) bnd2_->project(*(xs.get_2()));
107 }
108
120 Vector_SimOpt<Real> &xs = dynamic_cast<Vector_SimOpt<Real>&>(x);
121 if (bnd1_->isActivated()) bnd1_->projectInterior(*(xs.get_1()));
122 if (bnd2_->isActivated()) bnd2_->projectInterior(*(xs.get_2()));
123 }
124
136 void pruneUpperActive( Vector<Real> &v, const Vector<Real> &x, Real eps = Real(0) ) {
137 Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(v);
138 const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
139 if (bnd1_->isActivated()) bnd1_->pruneUpperActive(*(vs.get_1()),*(xs.get_1()),eps);
140 if (bnd2_->isActivated()) bnd2_->pruneUpperActive(*(vs.get_2()),*(xs.get_2()),eps);
141 }
142
156 void pruneUpperActive( Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real xeps = Real(0), Real geps = Real(0) ) {
157 Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(v);
158 const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(g);
159 const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
160 if (bnd1_->isActivated()) bnd1_->pruneUpperActive(*(vs.get_1()),*(gs.get_1()),*(xs.get_1()),xeps,geps);
161 if (bnd2_->isActivated()) bnd2_->pruneUpperActive(*(vs.get_2()),*(gs.get_2()),*(xs.get_2()),xeps,geps);
162 }
163
175 void pruneLowerActive( Vector<Real> &v, const Vector<Real> &x, Real eps = Real(0) ) {
176 Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(v);
177 const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
178 if (bnd1_->isActivated()) bnd1_->pruneLowerActive(*(vs.get_1()),*(xs.get_1()),eps);
179 if (bnd2_->isActivated()) bnd2_->pruneLowerActive(*(vs.get_2()),*(xs.get_2()),eps);
180 }
181
195 void pruneLowerActive( Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real xeps = Real(0), Real geps = Real(0) ) {
196 Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(v);
197 const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(g);
198 const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
199 if (bnd1_->isActivated()) bnd1_->pruneLowerActive(*(vs.get_1()),*(gs.get_1()),*(xs.get_1()),xeps,geps);
200 if (bnd2_->isActivated()) bnd2_->pruneLowerActive(*(vs.get_2()),*(gs.get_2()),*(xs.get_2()),xeps,geps);
201 }
202
203 const Ptr<const Vector<Real>> getLowerBound( void ) const {
204 const Ptr<const Vector<Real>> l1 = bnd1_->getLowerBound();
205 const Ptr<const Vector<Real>> l2 = bnd2_->getLowerBound();
206 return makePtr<Vector_SimOpt<Real>>( constPtrCast<Vector<Real>>(l1),
207 constPtrCast<Vector<Real>>(l2) );
208 }
209
210 const Ptr<const Vector<Real>> getUpperBound(void) const {
211 const Ptr<const Vector<Real>> u1 = bnd1_->getUpperBound();
212 const Ptr<const Vector<Real>> u2 = bnd2_->getUpperBound();
213 return makePtr<Vector_SimOpt<Real>>( constPtrCast<Vector<Real>>(u1),
214 constPtrCast<Vector<Real>>(u2) );
215 }
216
228 void pruneActive( Vector<Real> &v, const Vector<Real> &x, Real eps = Real(0) ) {
229 Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(v);
230 const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
231 if (bnd1_->isActivated()) bnd1_->pruneActive(*(vs.get_1()),*(xs.get_1()),eps);
232 if (bnd2_->isActivated()) bnd2_->pruneActive(*(vs.get_2()),*(xs.get_2()),eps);
233 }
234
247 void pruneActive( Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real xeps = Real(0), Real geps = Real(0) ) {
248 Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(v);
249 const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(g);
250 const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
251 if (bnd1_->isActivated()) bnd1_->pruneActive(*(vs.get_1()),*(gs.get_1()),*(xs.get_1()),xeps,geps);
252 if (bnd2_->isActivated()) bnd2_->pruneActive(*(vs.get_2()),*(gs.get_2()),*(xs.get_2()),xeps,geps);
253 }
254
260 bool isFeasible( const Vector<Real> &v ) {
261 const Vector_SimOpt<Real> &vs = dynamic_cast<const Vector_SimOpt<Real>&>(v);
262 return (bnd1_->isFeasible(*(vs.get_1()))) && (bnd2_->isFeasible(*(vs.get_2())));
263 }
264
279 Vector_SimOpt<Real> &dvs = dynamic_cast<Vector_SimOpt<Real>&>(dv);
280 const Vector_SimOpt<Real> &vs = dynamic_cast<const Vector_SimOpt<Real>&>(v);
281 const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
282 const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(g);
283 if (bnd1_->isActivated()) bnd1_->applyInverseScalingFunction(*(dvs.get_1()),*(vs.get_1()),*(xs.get_1()),*(gs.get_1()));
284 if (bnd2_->isActivated()) bnd2_->applyInverseScalingFunction(*(dvs.get_2()),*(vs.get_2()),*(xs.get_2()),*(gs.get_2()));
285 }
286
301 Vector_SimOpt<Real> &dvs = dynamic_cast<Vector_SimOpt<Real>&>(dv);
302 const Vector_SimOpt<Real> &vs = dynamic_cast<const Vector_SimOpt<Real>&>(v);
303 const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
304 const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(g);
305 if (bnd1_->isActivated()) bnd1_->applyScalingFunctionJacobian(*(dvs.get_1()),*(vs.get_1()),*(xs.get_1()),*(gs.get_1()));
306 if (bnd2_->isActivated()) bnd2_->applyScalingFunctionJacobian(*(dvs.get_2()),*(vs.get_2()),*(xs.get_2()),*(gs.get_2()));
307 }
308
309}; // class BoundConstraint_SimOpt
310
311} // namespace ROL
312
313#endif
Contains definitions of custom data types in ROL.
const Ptr< const Vector< Real > > getLowerBound(void) const
Return the ref count pointer to the lower bound vector.
bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
BoundConstraint_SimOpt(const Ptr< BoundConstraint< Real > > &bnd, bool optBnd=true)
Constructor for single bound constraint.
void pruneUpperActive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the upper -active set.
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the -active set.
BoundConstraint_SimOpt(const Ptr< BoundConstraint< Real > > &bnd, const Vector< Real > &x, bool optBnd=true)
Constructor for single bound constraint.
const Ptr< BoundConstraint< Real > > bnd1_
void project(Vector< Real > &x)
Project optimization variables onto the bounds.
void pruneLowerActive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the lower -active set.
void applyInverseScalingFunction(Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const
Apply inverse scaling function.
void projectInterior(Vector< Real > &x)
Project optimization variables into the interior of the feasible set.
const Ptr< BoundConstraint< Real > > bnd2_
const Ptr< const Vector< Real > > getUpperBound(void) const
Return the ref count pointer to the upper bound vector.
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 -binding set.
void applyScalingFunctionJacobian(Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const
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))
Set variables to zero if they correspond to the lower -binding set.
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 -binding set.
BoundConstraint_SimOpt(const Ptr< BoundConstraint< Real > > &bnd1, const Ptr< BoundConstraint< Real > > &bnd2)
Default constructor.
void deactivate(void)
Turn off bounds.
void activate(void)
Turn on bounds.
Defines the linear algebra or vector space interface for simulation-based optimization.
ROL::Ptr< const Vector< Real > > get_2() const
ROL::Ptr< const Vector< Real > > get_1() const
Defines the linear algebra or vector space interface.