ROL
ROL_BoundConstraint_SimOpt.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_BOUND_CONSTRAINT_SIMOPT_H
45 #define ROL_BOUND_CONSTRAINT_SIMOPT_H
46 
47 #include "ROL_BoundConstraint.hpp"
48 #include "ROL_Vector_SimOpt.hpp"
49 #include "ROL_Types.hpp"
50 #include <iostream>
51 
70 namespace ROL {
71 
72 template <class Real>
74 private:
75  Ptr<BoundConstraint<Real>> bnd1_;
76  Ptr<BoundConstraint<Real>> bnd2_;
77 
78 public:
80 
86  const Ptr<BoundConstraint<Real>> &bnd2)
87  : bnd1_(bnd1), bnd2_(bnd2) {
88  if ( bnd1_->isActivated() || bnd2_->isActivated() ) {
90  }
91  else {
93  }
94  }
95 
103  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
104  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(
105  dynamic_cast<const Vector<Real>&>(x));
106  if ( bnd1_->isActivated() ) {
107  bnd1_->update(*(xs.get_1()),flag,iter);
108  }
109  if ( bnd2_->isActivated() ) {
110  bnd2_->update(*(xs.get_2()),flag,iter);
111  }
112  }
113 
122  void project( Vector<Real> &x ) {
123  Vector_SimOpt<Real> &xs = dynamic_cast<Vector_SimOpt<Real>&>(
124  dynamic_cast<Vector<Real>&>(x));
125  if ( bnd1_->isActivated() ) {
126  bnd1_->project(*(xs.get_1()));
127  }
128  if ( bnd2_->isActivated() ) {
129  bnd2_->project(*(xs.get_2()));
130  }
131  }
132 
144  Vector_SimOpt<Real> &xs = dynamic_cast<Vector_SimOpt<Real>&>(
145  dynamic_cast<Vector<Real>&>(x));
146  if ( bnd1_->isActivated() ) {
147  bnd1_->projectInterior(*(xs.get_1()));
148  }
149  if ( bnd2_->isActivated() ) {
150  bnd2_->projectInterior(*(xs.get_2()));
151  }
152  }
153 
160  bool checkMultipliers( const Vector<Real> &l, const Vector<Real> &x ) {
161  const Vector_SimOpt<Real> &ls = dynamic_cast<const Vector_SimOpt<Real>&>(
162  dynamic_cast<const Vector<Real>&>(l));
163  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(
164  dynamic_cast<const Vector<Real>&>(x));
165  bool nn1 = true;
166  if ( bnd1_->isActivated() ) {
167  nn1 = bnd1_->checkMultipliers(*(ls.get_1()),*(xs.get_1()));
168  }
169  bool nn2 = true;
170  if ( bnd2_->isActivated() ) {
171  nn2 = bnd2_->checkMultipliers(*(ls.get_2()),*(xs.get_2()));
172  }
173  return (nn1 && nn2);
174  }
175 
187  void pruneUpperActive( Vector<Real> &v, const Vector<Real> &x, Real eps = Real(0) ) {
188  Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(
189  dynamic_cast<Vector<Real>&>(v));
190  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(
191  dynamic_cast<const Vector<Real>&>(x));
192  if ( bnd1_->isActivated() ) {
193  bnd1_->pruneUpperActive(*(vs.get_1()),*(xs.get_1()),eps);
194  }
195  if ( bnd2_->isActivated() ) {
196  bnd2_->pruneUpperActive(*(vs.get_2()),*(xs.get_2()),eps);
197  }
198  }
199 
213  void pruneUpperActive( Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real xeps = Real(0), Real geps = Real(0) ) {
214  Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(
215  dynamic_cast<Vector<Real>&>(v));
216  const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(
217  dynamic_cast<const Vector<Real>&>(g));
218  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(
219  dynamic_cast<const Vector<Real>&>(x));
220  if ( bnd1_->isActivated() ) {
221  bnd1_->pruneUpperActive(*(vs.get_1()),*(gs.get_1()),*(xs.get_1()),xeps,geps);
222  }
223  if ( bnd2_->isActivated() ) {
224  bnd2_->pruneUpperActive(*(vs.get_2()),*(gs.get_2()),*(xs.get_2()),xeps,geps);
225  }
226  }
227 
239  void pruneLowerActive( Vector<Real> &v, const Vector<Real> &x, Real eps = Real(0) ) {
240  Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(
241  dynamic_cast<Vector<Real>&>(v));
242  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(
243  dynamic_cast<const Vector<Real>&>(x));
244  if ( bnd1_->isActivated() ) {
245  bnd1_->pruneLowerActive(*(vs.get_1()),*(xs.get_1()),eps);
246  }
247  if ( bnd2_->isActivated() ) {
248  bnd2_->pruneLowerActive(*(vs.get_2()),*(xs.get_2()),eps);
249  }
250  }
251 
265  void pruneLowerActive( Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real xeps = Real(0), Real geps = Real(0) ) {
266  Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(
267  dynamic_cast<Vector<Real>&>(v));
268  const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(
269  dynamic_cast<const Vector<Real>&>(g));
270  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(
271  dynamic_cast<const Vector<Real>&>(x));
272  if ( bnd1_->isActivated() ) {
273  bnd1_->pruneLowerActive(*(vs.get_1()),*(gs.get_1()),*(xs.get_1()),xeps,geps);
274  }
275  if ( bnd2_->isActivated() ) {
276  bnd2_->pruneLowerActive(*(vs.get_2()),*(gs.get_2()),*(xs.get_2()),xeps,geps);
277  }
278  }
279 
280  const Ptr<const Vector<Real>> getLowerBound( void ) const {
281  const Ptr<const Vector<Real>> l1 = bnd1_->getLowerBound();
282  const Ptr<const Vector<Real>> l2 = bnd2_->getLowerBound();
283  return makePtr<Vector_SimOpt<Real>>( constPtrCast<Vector<Real>>(l1),
284  constPtrCast<Vector<Real>>(l2) );
285  }
286 
287  const Ptr<const Vector<Real>> getUpperBound(void) const {
288  const Ptr<const Vector<Real>> u1 = bnd1_->getUpperBound();
289  const Ptr<const Vector<Real>> u2 = bnd2_->getUpperBound();
290  return makePtr<Vector_SimOpt<Real>>( constPtrCast<Vector<Real>>(u1),
291  constPtrCast<Vector<Real>>(u2) );
292  }
293 
305  void pruneActive( Vector<Real> &v, const Vector<Real> &x, Real eps = Real(0) ) {
306  Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(
307  dynamic_cast<Vector<Real>&>(v));
308  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(
309  dynamic_cast<const Vector<Real>&>(x));
310  if ( bnd1_->isActivated() ) {
311  bnd1_->pruneActive(*(vs.get_1()),*(xs.get_1()),eps);
312  }
313  if ( bnd2_->isActivated() ) {
314  bnd2_->pruneActive(*(vs.get_2()),*(xs.get_2()),eps);
315  }
316  }
317 
330  void pruneActive( Vector<Real> &v, const Vector<Real> &g, const Vector<Real> &x, Real xeps = Real(0), Real geps = Real(0) ) {
331  Vector_SimOpt<Real> &vs = dynamic_cast<Vector_SimOpt<Real>&>(v);
332  const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(g);
333  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
334  if ( bnd1_->isActivated() ) {
335  bnd1_->pruneActive(*(vs.get_1()),*(gs.get_1()),*(xs.get_1()),xeps,geps);
336  }
337  if ( bnd2_->isActivated() ) {
338  bnd2_->pruneActive(*(vs.get_2()),*(gs.get_2()),*(xs.get_2()),xeps,geps);
339  }
340  }
341 
347  bool isFeasible( const Vector<Real> &v ) {
348  const Vector_SimOpt<Real> &vs = dynamic_cast<const Vector_SimOpt<Real>&>(v);
349  return (bnd1_->isFeasible(*(vs.get_1()))) && (bnd2_->isFeasible(*(vs.get_2())));
350  }
351 
365  void applyInverseScalingFunction(Vector<Real> &dv, const Vector<Real> &v, const Vector<Real> &x, const Vector<Real> &g) const{
366  Vector_SimOpt<Real> &dvs = dynamic_cast<Vector_SimOpt<Real>&>(dv);
367  const Vector_SimOpt<Real> &vs = dynamic_cast<const Vector_SimOpt<Real>&>(v);
368  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
369  const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(g);
370  if ( bnd1_->isActivated() ) {
371  bnd1_->applyInverseScalingFunction(*(dvs.get_1()),*(vs.get_1()),*(xs.get_1()),*(gs.get_1()));
372  }
373  if ( bnd2_->isActivated() ) {
374  bnd2_->applyInverseScalingFunction(*(dvs.get_2()),*(vs.get_2()),*(xs.get_2()),*(gs.get_2()));
375  }
376  }
377 
391  void applyScalingFunctionJacobian(Vector<Real> &dv, const Vector<Real> &v, const Vector<Real> &x, const Vector<Real> &g) const {
392  Vector_SimOpt<Real> &dvs = dynamic_cast<Vector_SimOpt<Real>&>(dv);
393  const Vector_SimOpt<Real> &vs = dynamic_cast<const Vector_SimOpt<Real>&>(v);
394  const Vector_SimOpt<Real> &xs = dynamic_cast<const Vector_SimOpt<Real>&>(x);
395  const Vector_SimOpt<Real> &gs = dynamic_cast<const Vector_SimOpt<Real>&>(g);
396  if ( bnd1_->isActivated() ) {
397  bnd1_->applyScalingFunctionJacobian(*(dvs.get_1()),*(vs.get_1()),*(xs.get_1()),*(gs.get_1()));
398  }
399  if ( bnd2_->isActivated() ) {
400  bnd2_->applyScalingFunctionJacobian(*(dvs.get_2()),*(vs.get_2()),*(xs.get_2()),*(gs.get_2()));
401  }
402  }
403 
404 }; // class BoundConstraint
405 
406 } // namespace ROL
407 
408 #endif
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the -active set.
bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
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 pruneLowerActive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the lower -active set.
Defines the linear algebra or vector space interface for simulation-based optimization.
void activate(void)
Turn on bounds.
Contains definitions of custom data types in ROL.
BoundConstraint_SimOpt(const Ptr< BoundConstraint< Real >> &bnd1, const Ptr< BoundConstraint< Real >> &bnd2)
Default constructor.
void applyScalingFunctionJacobian(Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const
Apply scaling function Jacobian.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
ROL::Ptr< const Vector< Real > > get_2() const
Ptr< BoundConstraint< Real > > bnd1_
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.
bool checkMultipliers(const Vector< Real > &l, const Vector< Real > &x)
Determine if a vector of Lagrange multipliers is nonnegative components.
ROL::Ptr< const Vector< Real > > get_1() const
void projectInterior(Vector< Real > &x)
Project optimization variables into the interior of the feasible set.
const Ptr< const Vector< Real > > getLowerBound(void) const
Return the ref count pointer to the lower bound vector.
void applyInverseScalingFunction(Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const
Apply inverse scaling function.
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 update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update bounds.
Provides the interface to apply upper and lower bound constraints.
Ptr< BoundConstraint< Real > > bnd2_
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 deactivate(void)
Turn off bounds.
void project(Vector< Real > &x)
Project optimization variables onto the bounds.