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