ROL
ROL_DoubleDogLeg_U.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_U_H
45 #define ROL_DOUBLEDOGLEG_U_H
46 
51 #include "ROL_TrustRegion_U.hpp"
52 #include "ROL_Types.hpp"
53 
54 namespace ROL {
55 
56 template<class Real>
57 class DoubleDogLeg_U : public TrustRegion_U<Real> {
58 private:
59 
60  Ptr<Vector<Real>> primal_, dual_;
61 
62 public:
63 
65 
66  void initialize(const Vector<Real> &x, const Vector<Real> &g) {
67  primal_ = x.clone();
68  dual_ = g.clone();
69  }
70 
71  void solve( Vector<Real> &s,
72  Real &snorm,
73  Real &pRed,
74  int &iflag,
75  int &iter,
76  const Real del,
77  TrustRegionModel_U<Real> &model ) {
78  Real tol = std::sqrt(ROL_EPSILON<Real>());
79  const Real one(1), zero(0), half(0.5), p2(0.2), p8(0.8), two(2);
80  // Set s to be the (projected) gradient
81  s.set(model.getGradient()->dual());
82  // Compute (quasi-)Newton step
83  model.invHessVec(*primal_,*model.getGradient(),s,tol);
84  Real sNnorm = primal_->norm();
85  Real tmp = -primal_->dot(s);
86  Real gsN = std::abs(tmp);
87  // Check if (quasi-)Newton step is feasible
88  if ( tmp >= zero ) {
89  // Use the Cauchy point
90  model.hessVec(*dual_,s,s,tol);
91  //Real gBg = dual_->dot(s.dual());
92  Real gBg = dual_->apply(s);
93  Real gnorm = s.dual().norm();
94  Real gg = gnorm*gnorm;
95  Real alpha = del/gnorm;
96  if ( gBg > ROL_EPSILON<Real>() ) {
97  alpha = std::min(gg/gBg, del/gnorm);
98  }
99  s.scale(-alpha);
100  snorm = alpha*gnorm;
101  iflag = 2;
102  pRed = alpha*(gg - half*alpha*gBg);
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(*primal_);
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(*dual_,s,s,tol);
115  Real alpha = zero;
116  Real beta = zero;
117  Real gnorm = s.norm();
118  Real gnorm2 = gnorm*gnorm;
119  //Real gBg = dual_->dot(s.dual());
120  Real gBg = dual_->apply(s);
121  Real gamma1 = gnorm/gBg;
122  Real gamma2 = gnorm/gsN;
123  Real eta = p8*gamma1*gamma2 + p2;
124  if (eta*sNnorm <= del || gBg <= zero) { // Dogleg Point is inside trust region
125  alpha = del/sNnorm;
126  beta = zero;
127  s.set(*primal_);
128  s.scale(-alpha);
129  snorm = del;
130  iflag = 1;
131  }
132  else {
133  if (gnorm2*gamma1 >= del) { // Cauchy Point is outside trust region
134  alpha = zero;
135  beta = -del/gnorm;
136  s.scale(beta);
137  snorm = del;
138  iflag = 2;
139  }
140  else { // Find convex combination of Cauchy and Dogleg point
141  s.scale(-gamma1*gnorm);
142  primal_->scale(eta);
143  primal_->plus(s);
144  primal_->scale(-one);
145  Real wNorm = primal_->dot(*primal_);
146  Real sigma = del*del-std::pow(gamma1*gnorm,two);
147  Real phi = s.dot(*primal_);
148  Real theta = (-phi + std::sqrt(phi*phi+wNorm*sigma))/wNorm;
149  s.axpy(theta,*primal_);
150  snorm = del;
151  alpha = theta*eta;
152  beta = (one-theta)*(-gamma1*gnorm);
153  iflag = 3;
154  }
155  }
156  pRed = -(alpha*(half*alpha-one)*gsN + half*beta*beta*gBg + beta*(one-alpha)*gnorm2);
157  }
158  }
159  }
160 };
161 
162 } // namespace ROL
163 
164 #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 invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply inverse Hessian approximation to vector.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:153
Contains definitions of custom data types in ROL.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
virtual Real dot(const Vector &x) const =0
Compute where .
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply Hessian approximation to vector.
void solve(Vector< Real > &s, Real &snorm, Real &pRed, int &iflag, int &iter, const Real del, TrustRegionModel_U< Real > &model)
Provides the interface to evaluate trust-region model functions.
virtual const Ptr< const Vector< Real > > getGradient(void) const
Ptr< Vector< Real > > dual_
Ptr< Vector< Real > > primal_
Provides interface for and implements trust-region subproblem solvers.
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
virtual Real norm() const =0
Returns where .
Provides interface for the double dog leg trust-region subproblem solver.
void initialize(const Vector< Real > &x, const Vector< Real > &g)