ROL
ROL_DoubleDogLeg.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_DOUBLEDOGLEG_H
45 #define ROL_DOUBLEDOGLEG_H
46 
51 #include "ROL_TrustRegion.hpp"
52 #include "ROL_Types.hpp"
53 
54 namespace ROL {
55 
56 template<class Real>
57 class DoubleDogLeg : public TrustRegion<Real> {
58 private:
59 
60  ROL::Ptr<CauchyPoint<Real> > cpt_;
61 
62  ROL::Ptr<Vector<Real> > s_;
63  ROL::Ptr<Vector<Real> > v_;
64  ROL::Ptr<Vector<Real> > Hp_;
65 
66  Real pRed_;
67 
68 public:
69 
70  // Constructor
71  DoubleDogLeg( ROL::ParameterList &parlist ) : TrustRegion<Real>(parlist), pRed_(0) {
72  cpt_ = ROL::makePtr<CauchyPoint<Real>>(parlist);
73  }
74 
75  void initialize( const Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g) {
77  cpt_->initialize(x,s,g);
78  s_ = s.clone();
79  v_ = s.clone();
80  Hp_ = g.clone();
81  }
82 
83  void run( Vector<Real> &s,
84  Real &snorm,
85  int &iflag,
86  int &iter,
87  const Real del,
88  TrustRegionModel<Real> &model ) {
89  Real tol = std::sqrt(ROL_EPSILON<Real>());
90  const Real one(1), zero(0), half(0.5), p2(0.2), p8(0.8), two(2);
91  // Set s to be the (projected) gradient
92  model.dualTransform(*Hp_,*model.getGradient());
93  s.set(Hp_->dual());
94  // Compute (quasi-)Newton step
95  model.invHessVec(*s_,*Hp_,s,tol);
96  Real sNnorm = s_->norm();
97  Real tmp = -s_->dot(s);
98  bool negCurv = (tmp > zero ? true : false);
99  Real gsN = std::abs(tmp);
100  // Check if (quasi-)Newton step is feasible
101  if ( negCurv ) {
102  // Use Cauchy point
103  cpt_->run(s,snorm,iflag,iter,del,model);
104  pRed_ = cpt_->getPredictedReduction();
105  iflag = 2;
106  }
107  else {
108  // Approximately solve trust region subproblem using double dogleg curve
109  if (sNnorm <= del) { // Use the (quasi-)Newton step
110  s.set(*s_);
111  s.scale(-one);
112  snorm = sNnorm;
113  pRed_ = half*gsN;
114  iflag = 0;
115  }
116  else { // The (quasi-)Newton step is outside of trust region
117  model.hessVec(*Hp_,s,s,tol);
118  Real alpha = zero;
119  Real beta = zero;
120  Real gnorm = s.norm();
121  Real gnorm2 = gnorm*gnorm;
122  Real gBg = Hp_->dot(s.dual());
123  Real gamma1 = gnorm/gBg;
124  Real gamma2 = gnorm/gsN;
125  Real eta = p8*gamma1*gamma2 + p2;
126  if (eta*sNnorm <= del || gBg <= zero) { // Dogleg Point is inside trust region
127  alpha = del/sNnorm;
128  beta = zero;
129  s.set(*s_);
130  s.scale(-alpha);
131  snorm = del;
132  iflag = 1;
133  }
134  else {
135  if (gnorm2*gamma1 >= del) { // Cauchy Point is outside trust region
136  alpha = zero;
137  beta = -del/gnorm;
138  s.scale(beta);
139  snorm = del;
140  iflag = 2;
141  }
142  else { // Find convex combination of Cauchy and Dogleg point
143  s.scale(-gamma1*gnorm);
144  v_->set(s);
145  v_->axpy(eta,*s_);
146  v_->scale(-one);
147  Real wNorm = v_->dot(*v_);
148  Real sigma = del*del-std::pow(gamma1*gnorm,two);
149  Real phi = s.dot(*v_);
150  Real theta = (-phi + std::sqrt(phi*phi+wNorm*sigma))/wNorm;
151  s.axpy(theta,*v_);
152  snorm = del;
153  alpha = theta*eta;
154  beta = (one-theta)*(-gamma1*gnorm);
155  iflag = 3;
156  }
157  }
158  pRed_ = -(alpha*(half*alpha-one)*gsN + half*beta*beta*gBg + beta*(one-alpha)*gnorm2);
159  }
160  }
161  model.primalTransform(*s_,s);
162  s.set(*s_);
163  snorm = s.norm();
165  }
166 };
167 
168 }
169 
170 #endif
ROL::Ptr< Vector< Real > > s_
void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
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
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
virtual Real dot(const Vector &x) const =0
Compute where .
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
Provides interface for the double dog leg trust-region subproblem solver.
void setPredictedReduction(const Real pRed)
ROL::Ptr< Vector< Real > > v_
DoubleDogLeg(ROL::ParameterList &parlist)
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.
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
virtual Real norm() const =0
Returns where .
void run(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)
ROL::Ptr< Vector< Real > > Hp_
ROL::Ptr< CauchyPoint< Real > > cpt_
virtual void primalTransform(Vector< Real > &tv, const Vector< Real > &v)