ROL
ROL_ColemanLiModel.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_COLEMANLIMODEL_HPP
45 #define ROL_COLEMANLIMODEL_HPP
46 
47 #include "ROL_TrustRegionModel.hpp"
48 
57 namespace ROL {
58 
59 template<class Real>
60 class ColemanLiModel : public TrustRegionModel<Real> {
61 private:
62  Ptr<Vector<Real>> prim_, dual_, hv_; // Auxiliary storage
63  Ptr<Vector<Real>> step_; // Step storage
64  Ptr<Vector<Real>> cauchyStep_, cauchyScal_; // Cauchy point vectors
65  Ptr<Vector<Real>> reflectStep_, reflectScal_; // Reflective step vectors
66  Ptr<Vector<Real>> Dmat_; // sqrt(abs(v))
67  Ptr<Vector<Real>> Cmat_; // diag(g) * dv/dx
68  Ptr<Vector<Real>> lx_, ux_; // Temporary vectors for bound computation
69 
70  Real TRradius_; // Trust-region radius
71  const Real stepBackMax_, stepBackScale_; // Primal transform parameters
72  const bool singleReflect_; // Use single reflection
73  Real sCs_, pred_; // Actual/predicted reduction
74 
75  Elementwise::Multiply<Real> mult_; // Elementwise multiply
76  Elementwise::Divide<Real> div_; // Elementwise division
77 
78  // Apply diagonal D matrix
79  void applyD( Vector<Real> &Dv, const Vector<Real> &v ) {
80  Dv.set(v);
81  Dv.applyBinary(div_,*Dmat_);
82  }
83 
84  // Apply inverse of diagonal D matrix
85  void applyInverseD( Vector<Real> &Dv, const Vector<Real> &v ) {
86  Dv.set(v);
87  Dv.applyBinary(mult_,*Dmat_);
88  }
89 
90  // Apply diagonal C matrix
91  void applyC( Vector<Real> &Cv, const Vector<Real> &v ) {
92  Cv.set(v);
93  Cv.applyBinary(mult_, *Cmat_);
94  }
95 
96  void constructC(void) {
97  const Ptr<const Vector<Real>> gc = TrustRegionModel<Real>::getGradient();
98  const Ptr<const Vector<Real>> l = TrustRegionModel<Real>::getBoundConstraint()->getLowerBound();
99  const Ptr<const Vector<Real>> u = TrustRegionModel<Real>::getBoundConstraint()->getUpperBound();
100 
101  // Set Cmat_ to be the sign of the gradient
102  Cmat_->set(gc->dual());
103  Cmat_->applyUnary(Elementwise::Sign<Real>());
104  // If g < 0 and u = INF then set Cmat_ to zero
105  class NegGradInfU : public Elementwise::BinaryFunction<Real> {
106  public:
107  NegGradInfU(void) {}
108  Real apply(const Real &x, const Real &y) const {
109  const Real zero(0), one(1), INF(ROL_INF<Real>());
110  return (x < zero && y == INF) ? zero : one;
111  }
112  };
113  prim_->set(gc->dual());
114  prim_->applyBinary(NegGradInfU(), *u);
115  Cmat_->applyBinary(mult_, *prim_);
116  // If g >= 0 and l = -INF then set Cmat_ to zero
117  class PosGradNinfL : public Elementwise::BinaryFunction<Real> {
118  public:
119  PosGradNinfL(void) {}
120  Real apply(const Real &x, const Real &y) const {
121  const Real zero(0), one(1), NINF(ROL_NINF<Real>());
122  return (x >= zero && y == NINF) ? zero : one;
123  }
124  };
125  prim_->set(gc->dual());
126  prim_->applyBinary(PosGradNinfL(), *l);
127  Cmat_->applyBinary(mult_, *prim_);
128  // Pointwise multiply Cmat_ with the gradient
129  Cmat_->applyBinary(mult_, gc->dual());
130  }
131 
132  void constructInverseD(void) {
133  const Ptr<const Vector<Real>> xc = TrustRegionModel<Real>::getIterate();
134  const Ptr<const Vector<Real>> gc = TrustRegionModel<Real>::getGradient();
135  const Ptr<const Vector<Real>> l = TrustRegionModel<Real>::getBoundConstraint()->getLowerBound();
136  const Ptr<const Vector<Real>> u = TrustRegionModel<Real>::getBoundConstraint()->getUpperBound();
137  const Real zero(0), one(1), INF(ROL_INF<Real>()), NINF(ROL_NINF<Real>());
138  const int LESS_THAN = 0;
139  const int EQUAL_TO = 1;
140  const int GREATER_THAN = 2;
141 
142  Dmat_->zero();
143  // CASE (i)
144  // Mask for negative gradient (m1 is 1 if g is negative and 0 otherwise)
145  reflectStep_->applyBinary(Elementwise::ValueSet<Real>(zero, LESS_THAN),gc->dual());
146  // Mask for finite upper bounds (m2 is 1 if u is finite and 0 otherwise)
147  reflectScal_->applyBinary(Elementwise::ValueSet<Real>(INF, LESS_THAN),*u);
148  // Mask for g_i < 0 and u_i < inf
149  reflectScal_->applyBinary(mult_,*reflectStep_);
150  // prim_i = { u_i-x_i if g_i < 0 and u_i < inf
151  // { 0 otherwise
152  prim_->set(*u); prim_->axpy(-one,*xc);
153  prim_->applyBinary(mult_,*reflectScal_);
154  // Add to D
155  Dmat_->plus(*prim_);
156 
157  // CASE (iii)
158  // Mask for infinite upper bounds
159  reflectScal_->applyBinary(Elementwise::ValueSet<Real>(INF, EQUAL_TO),*u);
160  // Mask for g_i < 0 and u_i = inf
161  reflectScal_->applyBinary(mult_,*reflectStep_);
162  // prim_i = { -1 if g_i < 0 and u_i = inf
163  // { 0 otherwise
164  prim_->applyUnary(Elementwise::Fill<Real>(-one));
165  prim_->applyBinary(mult_,*reflectScal_);
166  // Add to D
167  Dmat_->plus(*prim_);
168 
169  // CASE (ii)
170  // m1 = 1-m1
171  reflectStep_->scale(-one);
172  reflectStep_->applyUnary(Elementwise::Shift<Real>(one));
173  // Mask for finite lower bounds
174  reflectScal_->applyBinary(Elementwise::ValueSet<Real>(NINF, GREATER_THAN),*l);
175  // Zero out elements of Jacobian with l_i=-inf
176  reflectScal_->applyBinary(mult_,*reflectStep_);
177  // prim_i = { x_i-l_i if g_i >= 0 and l_i > -inf
178  // { 0 otherwise
179  prim_->set(*xc); prim_->axpy(-one,*l);
180  prim_->applyBinary(mult_,*reflectScal_);
181  // Add to D
182  Dmat_->plus(*prim_);
183 
184  // CASE (iv)
185  // Mask for infinite lower bounds
186  reflectScal_->applyBinary(Elementwise::ValueSet<Real>(NINF, EQUAL_TO),*l);
187  // Mask for g_i>=0 and l_i=-inf
188  reflectScal_->applyBinary(mult_,*reflectStep_);
189  // prim_i = { 1 if g_i >= 0 and l_i = -inf
190  // { 0 otherwise
191  prim_->applyUnary(Elementwise::Fill<Real>(one));
192  prim_->applyBinary(mult_,*reflectScal_);
193  // Add to D
194  Dmat_->plus(*prim_);
195 
196  // d_i = { u_i-x_i if g_i < 0, u_i<inf
197  // { -1 if g_i < 0, u_i=inf
198  // { x_i-l_i if g_i >= 0, l_i>-inf
199  // { 1 if g_i >= 0, l_i=-inf
200  Dmat_->applyUnary(Elementwise::AbsoluteValue<Real>());
201  Dmat_->applyUnary(Elementwise::SquareRoot<Real>());
202  }
203 
204  // Build diagonal D and C matrices
205  void initialize(const Vector<Real> &x, const Vector<Real> &g) {
206  prim_ = x.clone();
207  dual_ = g.clone();
208  hv_ = g.clone();
209  step_ = x.clone();
210  Dmat_ = x.clone();
211  Cmat_ = x.clone();
212  lx_ = x.clone();
213  ux_ = x.clone();
214 
215  cauchyStep_ = x.clone();
216  cauchyScal_ = x.clone();
217  reflectStep_ = x.clone();
218  reflectScal_ = x.clone();
219  }
220 
221  public:
222 
224  const Vector<Real> &x, const Vector<Real> &g,
225  const Real stepBackMax = 0.9999, const Real stepBackScale = 1.0,
226  const bool singleReflect = true,
227  const Ptr<Secant<Real>> &secant = nullPtr,
228  const bool useSecantPrecond = false, const bool useSecantHessVec = false)
229  : TrustRegionModel<Real>::TrustRegionModel(obj,bnd,x,g,secant,useSecantPrecond,useSecantHessVec),
230  TRradius_(1), stepBackMax_(stepBackMax), stepBackScale_(stepBackScale),
231  singleReflect_(singleReflect), sCs_(0), pred_(0) {
232  initialize(x,g);
233  }
234 
236  const Vector<Real> &x, const Vector<Real> &g,
237  const Ptr<Secant<Real>> &secant = nullPtr) {
238  TrustRegionModel<Real>::update(obj,bnd,x,g,secant);
239  constructC();
241  }
242 
243  void setRadius(const Real del) {
244  TRradius_ = del;
245  }
246 
247  /***************************************************************************/
248  /********* BEGIN OBJECTIVE FUNCTION DEFINITIONS **************************/
249  /***************************************************************************/
250  // Note that s is the \f$\hat{s}\f$ and \f$\psi\f$ is the $\hat\psi$ from the paper
251  Real value( const Vector<Real> &s, Real &tol ) {
252  const Ptr<const Vector<Real>> gc = TrustRegionModel<Real>::getGradient();
253  // Apply Hessian to s
254  hessVec(*hv_, s, s, tol);
255  hv_->scale(static_cast<Real>(0.5));
256  // Form inv(D) * g
257  applyInverseD(*prim_, gc->dual());
258  // Add scaled gradient to Hessian in direction s
259  hv_->plus(prim_->dual());
260  return hv_->dot(s.dual());
261  }
262 
263  void gradient( Vector<Real> &g, const Vector<Real> &s, Real &tol ) {
264  const Ptr<const Vector<Real>> gc = TrustRegionModel<Real>::getGradient();
265  hessVec(g, s, s, tol);
266  applyInverseD(*prim_, gc->dual());
267  g.plus(prim_->dual());
268  }
269 
270  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &s, Real &tol ) {
271  const Ptr<const Vector<Real>> gc = TrustRegionModel<Real>::getGradient();
272  // Build B = inv(D) * Hessian * inv(D)
273  applyInverseD(*prim_, v);
275  applyInverseD(hv, *dual_);
276  // Build C = diag(g) J
277  applyC(*prim_, v);
278  hv.plus(prim_->dual());
279  }
280  /***************************************************************************/
281  /********* END OBJECTIVE FUNCTION DEFINITIONS ****************************/
282  /***************************************************************************/
283 
284  void dualTransform( Vector<Real> &tv, const Vector<Real> &v ) {
285  applyInverseD(tv, v);
286  }
287 
288  void primalTransform( Vector<Real> &tiv, const Vector<Real> &v ) {
289  Real tol = std::sqrt(ROL_EPSILON<Real>());
290 
291  /**************************************************************************/
292  /* PERFORM OPTIMAL SCALING OF TRUST REGION SUBPROBLEM SOLUTION */
293  /**************************************************************************/
294  applyInverseD(tiv, v);
295  // Get bounds on scalar variable
296  Real lowerBoundV(ROL_NINF<Real>()), upperBoundV(ROL_INF<Real>());
297  getScalarBounds(lowerBoundV, upperBoundV, tiv);
298  // Minimize one dimensional quadratic over bounds
299  Real tauV(1);
300  Real valueV = minimize1D(tauV, lowerBoundV, upperBoundV, v);
301 
302  /**************************************************************************/
303  /* COMPUTE CAUCHY POINT: STORED IN cauchyStep_ AND cauchyScal_ */
304  /**************************************************************************/
305  Real valueG = computeCauchyPoint();
306 
307  /**************************************************************************/
308  /* COMPUTE REFLECTIVE STEP: STORED IN reflectStep_ AND reflectScal_ */
309  /**************************************************************************/
310  if ( singleReflect_ ) {
312  }
313  else {
315  }
317  // Get bounds on scalar variable
318  Real lowerBoundR(ROL_NINF<Real>()), upperBoundR(ROL_INF<Real>());
319  getScalarBounds(lowerBoundR, upperBoundR, *reflectScal_);
320  // Minimize one dimensional quadratic over bounds
321  Real tauR(1);
322  Real valueR = minimize1D(tauR, lowerBoundR, upperBoundR, *reflectStep_);
323 
324  /**************************************************************************/
325  /* CHOOSE STEP THAT GIVES MOST PREDICTED REDUCTION */
326  /**************************************************************************/
327  Real VALUE(0);
328  bool useCauchyPoint = (valueG < valueV);
329  if (useCauchyPoint) {
330  VALUE = valueG;
331  tiv.set(*cauchyScal_);
332  // Store unscaled step
333  step_->set(*cauchyStep_);
334  }
335  else {
336  VALUE = valueV;
337  tiv.scale(tauV);
338  // Store unscaled step
339  step_->set(v);
340  step_->scale(tauV);
341  }
342  bool useReflectStep = (valueR < VALUE);
343  if (useReflectStep) {
344  VALUE = valueR;
345  tiv.set(*reflectScal_);
346  tiv.scale(tauR);
347  // Store unscaled step
348  step_->set(*reflectStep_);
349  step_->scale(tauR);
350  }
351 
352  /**************************************************************************/
353  /* ENSURE CHOSEN STEP IS STRICTLY FEASIBLE */
354  /**************************************************************************/
355  // Computed predicted reduction based on input step
356  if ( !isStrictlyFeasibleStep(tiv) ) {
357  Real snorm = step_->norm();
358  Real theta = std::max( stepBackMax_, static_cast<Real>(1) - stepBackScale_ * snorm);
359  tiv.scale(theta);
360  step_->scale(theta);
361  // Compute predicted reduction
362  pred_ = -value(*step_,tol);
363  }
364  else { // Theta is one
365  // Compute predicted reduction
366  pred_ = -VALUE;
367  }
368 
369  // Compute update for actual reduction
370  applyC(*prim_, *step_);
371  sCs_ = static_cast<Real>(-0.5) * prim_->dot(*step_);
372  }
373 
374  void updatePredictedReduction(Real &pred, const Vector<Real> &s) {
375  pred = pred_;
376  }
377 
378  void updateActualReduction(Real &ared, const Vector<Real> &s) {
379  ared += sCs_;
380  }
381 
382 private:
383 
384  void getScalarBounds( Real &lowerBound, Real &upperBound, const Vector<Real> &p ) {
385  const Ptr<const Vector<Real>> xc = TrustRegionModel<Real>::getIterate();
386  const Ptr<const Vector<Real>> l = TrustRegionModel<Real>::getBoundConstraint()->getLowerBound();
387  const Ptr<const Vector<Real>> u = TrustRegionModel<Real>::getBoundConstraint()->getUpperBound();
388  const Real one(1);
389  Real pnorm = p.norm();
390 
391  // Define elementwise functions
392  class PruneNegative : public Elementwise::BinaryFunction<Real> {
393  private:
394  const Real val_;
395  public:
396  PruneNegative( const Real val ) : val_(val) {}
397  Real apply(const Real &x, const Real &y) const {
398  return (y < static_cast<Real>(0)) ? x/y : val_;
399  }
400  };
401  class PrunePositive : public Elementwise::BinaryFunction<Real> {
402  private:
403  const Real val_;
404  public:
405  PrunePositive( const Real val ) : val_(val) {}
406  Real apply(const Real &x, const Real &y) const {
407  return (y > static_cast<Real>(0)) ? x/y : val_;
408  }
409  };
410 
411  // Max of (l-x)/p if p > 0
412  prim_->set(*l); prim_->axpy(-one,*xc);
413  prim_->applyBinary(PrunePositive(ROL_NINF<Real>()),p);
414  Real lowerBound1 = prim_->reduce(Elementwise::ReductionMax<Real>());
415  // Max of (u-x)/p if p < 0
416  prim_->set(*u); prim_->axpy(-one,*xc);
417  prim_->applyBinary(PruneNegative(ROL_NINF<Real>()),p);
418  Real lowerBound2 = prim_->reduce(Elementwise::ReductionMax<Real>());
419  // Lower bound
420  Real lowerBound3 = std::max(lowerBound1, lowerBound2);
421 
422  // Min of (u-x)/p if p > 0
423  prim_->set(*u); prim_->axpy(-one,*xc);
424  prim_->applyBinary(PrunePositive(ROL_INF<Real>()),p);
425  Real upperBound1 = prim_->reduce(Elementwise::ReductionMin<Real>());
426  // Max of (l-x)/p if p < 0
427  prim_->set(*l); prim_->axpy(-one,*xc);
428  prim_->applyBinary(PruneNegative(ROL_INF<Real>()),p);
429  Real upperBound2 = prim_->reduce(Elementwise::ReductionMin<Real>());
430  // Upper bound
431  Real upperBound3 = std::min(upperBound1, upperBound2);
432 
433  // Adjust for trust-region constraint
434  lowerBound = std::max(lowerBound3, -TRradius_/pnorm);
435  upperBound = std::min(upperBound3, TRradius_/pnorm);
436  }
437 
438  Real minimize1D(Real &tau, const Real lowerBound, const Real upperBound, const Vector<Real> &p) {
439  const Ptr<const Vector<Real>> gc = TrustRegionModel<Real>::getGradient();
440  Real tol = std::sqrt(ROL_EPSILON<Real>());
441 
442  // Compute coefficients of one dimensional quadratic
443  hessVec(*hv_, p, p, tol);
444  Real c2 = static_cast<Real>(0.5) * hv_->dot(p.dual());
445  applyInverseD(*prim_, gc->dual());
446  Real c1 = prim_->dot(p);
447 
448  // Minimize one dimensional quadratic over bounds
449  Real lval = (c2 * lowerBound + c1) * lowerBound;
450  Real rval = (c2 * upperBound + c1) * upperBound;
451  tau = (lval < rval) ? lowerBound : upperBound;
452  if (c2 > static_cast<Real>(0)) {
453  Real uncMin = static_cast<Real>(-0.5) * c1/c2;
454  tau = (uncMin > lowerBound && uncMin < upperBound) ? uncMin : tau;
455  }
456 
457  // Return minimal function value
458  return (c2 * tau + c1) * tau;
459  }
460 
461  Real computeCauchyPoint(void) {
462  // Set step = -inv(D) g
463  const Ptr<const Vector<Real>> gc = TrustRegionModel<Real>::getGradient();
464  applyInverseD(*cauchyStep_, gc->dual());
465  cauchyStep_->scale(static_cast<Real>(-1));
466 
467  // Scale Cauchy point
469 
470  // Scalar bounds
471  Real lowerBound(ROL_NINF<Real>()), upperBound(ROL_INF<Real>());
472  getScalarBounds(lowerBound, upperBound, *cauchyScal_);
473 
474  // Minimize 1D quadratic over bounds
475  Real tau(1), value(0);
476  value = minimize1D(tau, lowerBound, upperBound, *cauchyStep_);
477 
478  // Scale Cauchy point and return minimal function value
479  cauchyStep_->scale(tau);
480  cauchyScal_->scale(tau);
481  return value;
482  }
483 
485  const Ptr<const Vector<Real>> xc = TrustRegionModel<Real>::getIterate();
486  Real alpha = computeAlpha(Dv);
487  Rv.set(v);
488 
489  class LowerBound : public Elementwise::BinaryFunction<Real> {
490  public:
491  Real apply( const Real &x, const Real &y ) const {
492  return (x == y) ? static_cast<Real>(-1) : static_cast<Real>(1);
493  }
494  };
495  prim_->set(*xc); prim_->axpy(alpha,Dv);
496  prim_->applyBinary(LowerBound(),*TrustRegionModel<Real>::getBoundConstraint()->getLowerBound());
497  Rv.applyBinary(mult_,*prim_);
498 
499  class UpperBound : public Elementwise::BinaryFunction<Real> {
500  public:
501  Real apply( const Real &x, const Real &y ) const {
502  return (x == y) ? static_cast<Real>(-1) : static_cast<Real>(1);
503  }
504  };
505  prim_->set(*xc); prim_->axpy(alpha,Dv);
506  prim_->applyBinary(UpperBound(),*TrustRegionModel<Real>::getBoundConstraint()->getUpperBound());
507  Rv.applyBinary(mult_,*prim_);
508  }
509 
511  const Ptr<const Vector<Real>> xc = TrustRegionModel<Real>::getIterate();
512  Rv.set(v);
513 
514  class LowerBound : public Elementwise::BinaryFunction<Real> {
515  public:
516  Real apply( const Real &x, const Real &y ) const {
517  return (x < y) ? static_cast<Real>(-1) : static_cast<Real>(1);
518  }
519  };
520  prim_->set(*xc); prim_->plus(Dv);
521  prim_->applyBinary(LowerBound(),*TrustRegionModel<Real>::getBoundConstraint()->getLowerBound());
522  Rv.applyBinary(mult_,*prim_);
523 
524  class UpperBound : public Elementwise::BinaryFunction<Real> {
525  public:
526  Real apply( const Real &x, const Real &y ) const {
527  return (x > y) ? static_cast<Real>(-1) : static_cast<Real>(1);
528  }
529  };
530  prim_->set(*xc); prim_->plus(Dv);
531  prim_->applyBinary(UpperBound(),*TrustRegionModel<Real>::getBoundConstraint()->getUpperBound());
532  Rv.applyBinary(mult_,*prim_);
533  }
534 
535  Real computeAlpha( const Vector<Real> &d ) {
536  const Ptr<const Vector<Real>> xc = TrustRegionModel<Real>::getIterate();
537  const Real one(1);
538 
539  // Define elementwise functions
540  class SafeDivide : public Elementwise::BinaryFunction<Real> {
541  private:
542  const Real val_;
543  public:
544  SafeDivide( const Real val ) : val_(val) {}
545  Real apply(const Real &x, const Real &y) const {
546  const Real zero(0);
547  return (y == zero) ? val_ : x/y;
548  }
549  };
550 
551  // (l - x) / d
552  lx_->set(*TrustRegionModel<Real>::getBoundConstraint()->getLowerBound());
553  lx_->axpy(-one, *xc);
554  lx_->applyBinary(SafeDivide(ROL_INF<Real>()), d);
555 
556  // (u - x) / d
557  ux_->set(*TrustRegionModel<Real>::getBoundConstraint()->getUpperBound());
558  ux_->axpy(-one, *xc);
559  ux_->applyBinary(SafeDivide(ROL_INF<Real>()), d);
560 
561  // max{ (l - x) / d, (u - x) / d }
562  lx_->applyBinary(Elementwise::Max<Real>(),*ux_);
563 
564  // min{ max{ (l - x) / d, (u - x) / d } }
565  return lx_->reduce(Elementwise::ReductionMin<Real>());
566  }
567 
568  bool isStrictlyFeasibleStep( const Vector<Real> &d ) const {
569  const Ptr<const Vector<Real>> xc = TrustRegionModel<Real>::getIterate();
570 
571  class Greater : public Elementwise::BinaryFunction<Real> {
572  public:
573  Greater() {}
574  Real apply(const Real &x, const Real &y) const {
575  return (x > y) ? static_cast<Real>(1) : static_cast<Real>(0);
576  }
577  };
578  prim_->set(*xc); prim_->plus(d);
579  prim_->applyBinary(Greater(),*TrustRegionModel<Real>::getBoundConstraint()->getLowerBound());
580  Real lowerFeasible = prim_->reduce(Elementwise::ReductionMin<Real>());
581 
582  class Lesser : public Elementwise::BinaryFunction<Real> {
583  public:
584  Lesser() {}
585  Real apply(const Real &x, const Real &y) const {
586  return (x < y) ? static_cast<Real>(1) : static_cast<Real>(0);
587  }
588  };
589  prim_->set(*xc); prim_->plus(d);
590  prim_->applyBinary(Lesser(),*TrustRegionModel<Real>::getBoundConstraint()->getUpperBound());
591  Real upperFeasible = prim_->reduce(Elementwise::ReductionMin<Real>());
592 
593  return (upperFeasible * lowerFeasible > 0);
594  }
595 
596 }; // class ColemanLiModel
597 
598 }
599 
600 #endif // ROL_COLEMANLIMODEL_HPP
Provides the interface to evaluate objective functions.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: ROL_Vector.hpp:226
Ptr< Vector< Real > > lx_
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void plus(const Vector &x)=0
Compute , where .
void computeFullReflectiveStep(Vector< Real > &Rv, const Vector< Real > &v, const Vector< Real > &Dv)
Ptr< Vector< Real > > dual_
Ptr< Vector< Real > > prim_
Ptr< Vector< Real > > cauchyScal_
void initialize(const Vector< Real > &x, const Vector< Real > &g)
virtual void applyBinary(const Elementwise::BinaryFunction< Real > &f, const Vector &x)
Definition: ROL_Vector.hpp:236
void gradient(Vector< Real > &g, const Vector< Real > &s, Real &tol)
Compute gradient.
void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
void setRadius(const Real del)
void updatePredictedReduction(Real &pred, const Vector< Real > &s)
void updateActualReduction(Real &ared, const Vector< Real > &s)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply Hessian approximation to vector.
virtual void update(Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real >> &secant=nullPtr)
void primalTransform(Vector< Real > &tiv, const Vector< Real > &v)
Ptr< Vector< Real > > Cmat_
Ptr< Vector< Real > > Dmat_
Provides the interface to evaluate trust-region model functions.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
virtual const Ptr< const Vector< Real > > getGradient(void) const
Elementwise::Divide< Real > div_
Elementwise::Multiply< Real > mult_
void update(Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Ptr< Secant< Real >> &secant=nullPtr)
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
void applyD(Vector< Real > &Dv, const Vector< Real > &v)
bool isStrictlyFeasibleStep(const Vector< Real > &d) const
ColemanLiModel(Objective< Real > &obj, BoundConstraint< Real > &bnd, const Vector< Real > &x, const Vector< Real > &g, const Real stepBackMax=0.9999, const Real stepBackScale=1.0, const bool singleReflect=true, const Ptr< Secant< Real >> &secant=nullPtr, const bool useSecantPrecond=false, const bool useSecantHessVec=false)
void applyHessian(Vector< Real > &hv, const Vector< Real > &v, Real &tol)
Ptr< Vector< Real > > cauchyStep_
Ptr< Vector< Real > > ux_
Real minimize1D(Real &tau, const Real lowerBound, const Real upperBound, const Vector< Real > &p)
Real value(const Vector< Real > &s, Real &tol)
Compute value.
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:70
void applyInverseD(Vector< Real > &Dv, const Vector< Real > &v)
void applyC(Vector< Real > &Cv, const Vector< Real > &v)
Ptr< Vector< Real > > hv_
Provides the interface to evaluate interior trust-region model functions from the Coleman-Li bound co...
void getScalarBounds(Real &lowerBound, Real &upperBound, const Vector< Real > &p)
Real computeAlpha(const Vector< Real > &d)
Provides the interface to apply upper and lower bound constraints.
Ptr< Vector< Real > > reflectScal_
ROL::DiagonalOperator apply
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
virtual Real norm() const =0
Returns where .
Ptr< Vector< Real > > step_
virtual const Ptr< BoundConstraint< Real > > getBoundConstraint(void) const
Ptr< Vector< Real > > reflectStep_
void computeReflectiveStep(Vector< Real > &Rv, const Vector< Real > &v, const Vector< Real > &Dv)
virtual const Ptr< const Vector< Real > > getIterate(void) const