ROL
ROL_DogLeg.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_DOGLEG_H
45 #define ROL_DOGLEG_H
46 
51 #include "ROL_TrustRegion.hpp"
52 #include "ROL_Types.hpp"
53 
54 namespace ROL {
55 
56 template<class Real>
57 class DogLeg : public TrustRegion<Real> {
58 private:
59 
60  ROL::Ptr<CauchyPoint<Real> > cpt_;
61 
62  ROL::Ptr<Vector<Real> > s_;
63  ROL::Ptr<Vector<Real> > Hp_;
64 
65  Real pRed_;
66 
67 public:
68 
69  // Constructor
70  DogLeg( ROL::ParameterList &parlist ) : TrustRegion<Real>(parlist), pRed_(0) {
71  cpt_ = ROL::makePtr<CauchyPoint<Real>>(parlist);
72  }
73 
74  void initialize( const Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g) {
76  cpt_->initialize(x,s,g);
77  s_ = s.clone();
78  Hp_ = g.clone();
79  }
80 
81  void run( Vector<Real> &s,
82  Real &snorm,
83  int &iflag,
84  int &iter,
85  const Real del,
86  TrustRegionModel<Real> &model ) {
87  Real tol = std::sqrt(ROL_EPSILON<Real>());
88  const Real zero(0), half(0.5), one(1), two(2);
89  // Set s to be the (projected) gradient
90  model.dualTransform(*Hp_,*model.getGradient());
91  s.set(Hp_->dual());
92  // Compute (quasi-)Newton step
93  model.invHessVec(*s_,*Hp_,s,tol);
94  Real sNnorm = s_->norm();
95  Real gsN = -s_->dot(s);
96  bool negCurv = (gsN > zero ? true : false);
97  // Check if (quasi-)Newton step is feasible
98  if ( negCurv ) {
99  // Use Cauchy point
100  cpt_->run(s,snorm,iflag,iter,del,model);
101  pRed_ = cpt_->getPredictedReduction();
102  iflag = 2;
103  }
104  else {
105  // Approximately solve trust region subproblem using double dogleg curve
106  if (sNnorm <= del) { // Use the (quasi-)Newton step
107  s.set(*s_);
108  s.scale(-one);
109  snorm = sNnorm;
110  pRed_ = -half*gsN;
111  iflag = 0;
112  }
113  else { // The (quasi-)Newton step is outside of trust region
114  model.hessVec(*Hp_,s,s,tol);
115  Real alpha = zero;
116  Real beta = zero;
117  Real gnorm = s.norm();
118  Real gnorm2 = gnorm*gnorm;
119  Real gBg = Hp_->dot(s.dual());
120  Real gamma = gnorm2/gBg;
121  if ( gamma*gnorm >= del || gBg <= zero ) {
122  // Use Cauchy point
123  alpha = zero;
124  beta = del/gnorm;
125  s.scale(-beta);
126  snorm = del;
127  iflag = 2;
128  }
129  else {
130  // Use a convex combination of Cauchy point and (quasi-)Newton step
131  Real a = sNnorm*sNnorm + two*gamma*gsN + gamma*gamma*gnorm2;
132  Real b = -gamma*gsN - gamma*gamma*gnorm2;
133  Real c = gamma*gamma*gnorm2 - del*del;
134  alpha = (-b + sqrt(b*b - a*c))/a;
135  beta = gamma*(one-alpha);
136  s.scale(-beta);
137  s.axpy(-alpha,*s_);
138  snorm = del;
139  iflag = 1;
140  }
141  pRed_ = (alpha*(half*alpha-one)*gsN - half*beta*beta*gBg + beta*(one-alpha)*gnorm2);
142  }
143  }
144  model.primalTransform(*s_,s);
145  s.set(*s_);
146  snorm = s.norm();
147  // Update predicted reduction
149  }
150 };
151 
152 }
153 
154 #endif
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
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 axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:153
ROL::Ptr< Vector< Real > > Hp_
Definition: ROL_DogLeg.hpp:63
virtual void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
Contains definitions of custom data types in ROL.
Provides interface for and implements trust-region subproblem solvers.
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
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
void run(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)
Definition: ROL_DogLeg.hpp:81
void setPredictedReduction(const Real pRed)
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply Hessian approximation to vector.
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol)
Apply inverse Hessian approximation to vector.
ROL::Ptr< CauchyPoint< Real > > cpt_
Definition: ROL_DogLeg.hpp:60
ROL::Ptr< Vector< Real > > s_
Definition: ROL_DogLeg.hpp:62
Provides interface for dog leg trust-region subproblem solver.
Definition: ROL_DogLeg.hpp:57
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
virtual Real norm() const =0
Returns where .
DogLeg(ROL::ParameterList &parlist)
Definition: ROL_DogLeg.hpp:70
void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
Definition: ROL_DogLeg.hpp:74
virtual void primalTransform(Vector< Real > &tv, const Vector< Real > &v)