ROL
ROL_LineSearch.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_LINESEARCH_H
45 #define ROL_LINESEARCH_H
46 
51 #include "ROL_Ptr.hpp"
52 #include "ROL_ParameterList.hpp"
53 #include "ROL_Types.hpp"
54 #include "ROL_Vector.hpp"
55 #include "ROL_Objective.hpp"
56 #include "ROL_BoundConstraint.hpp"
57 #include "ROL_ScalarFunction.hpp"
58 
59 namespace ROL {
60 
61 template<class Real>
62 class LineSearch {
63 private:
64 
67 
68  bool useralpha_;
69  bool usePrevAlpha_; // Use the previous step's accepted alpha as an initial guess
70  Real alpha0_;
71  Real alpha0bnd_; // Lower bound for initial alpha...if below, set initial alpha to one
72  int maxit_;
73  Real c1_;
74  Real c2_;
75  Real c3_;
76  Real eps_;
77  Real fmin_; // smallest fval encountered
78  Real alphaMin_; // Alpha that yields the smallest fval encountered
79  bool acceptMin_; // Use smallest fval if sufficient decrease not satisfied
80  bool itcond_; // true if maximum function evaluations reached
81 
82  ROL::Ptr<Vector<Real> > xtst_;
83  ROL::Ptr<Vector<Real> > d_;
84  ROL::Ptr<Vector<Real> > g_;
85  ROL::Ptr<Vector<Real> > grad_;
86 // ROL::Ptr<const Vector<Real> > grad_;
87 
88 public:
89 
90 
91  virtual ~LineSearch() {}
92 
93  // Constructor
94  LineSearch( ROL::ParameterList &parlist ) : eps_(0) {
95  Real one(1), p9(0.9), p6(0.6), p4(0.4), oem4(1.e-4), zero(0);
96  // Enumerations
97  edesc_ = StringToEDescent(parlist.sublist("Step").sublist("Line Search").sublist("Descent Method").get("Type","Quasi-Newton Method"));
98  econd_ = StringToECurvatureCondition(parlist.sublist("Step").sublist("Line Search").sublist("Curvature Condition").get("Type","Strong Wolfe Conditions"));
99  // Linesearch Parameters
100  alpha0_ = parlist.sublist("Step").sublist("Line Search").get("Initial Step Size",one);
101  alpha0bnd_ = parlist.sublist("Step").sublist("Line Search").get("Lower Bound for Initial Step Size",one);
102  useralpha_ = parlist.sublist("Step").sublist("Line Search").get("User Defined Initial Step Size",false);
103  usePrevAlpha_ = parlist.sublist("Step").sublist("Line Search").get("Use Previous Step Length as Initial Guess",false);
104  acceptMin_ = parlist.sublist("Step").sublist("Line Search").get("Accept Linesearch Minimizer",false);
105  maxit_ = parlist.sublist("Step").sublist("Line Search").get("Function Evaluation Limit",20);
106  c1_ = parlist.sublist("Step").sublist("Line Search").get("Sufficient Decrease Tolerance",oem4);
107  c2_ = parlist.sublist("Step").sublist("Line Search").sublist("Curvature Condition").get("General Parameter",p9);
108  c3_ = parlist.sublist("Step").sublist("Line Search").sublist("Curvature Condition").get("Generalized Wolfe Parameter",p6);
109 
110  fmin_ = std::numeric_limits<Real>::max();
111  alphaMin_ = 0;
112  itcond_ = false;
113 
114  c1_ = ((c1_ < zero) ? oem4 : c1_);
115  c2_ = ((c2_ < zero) ? p9 : c2_);
116  c3_ = ((c3_ < zero) ? p9 : c3_);
117  if ( c2_ <= c1_ ) {
118  c1_ = oem4;
119  c2_ = p9;
120  }
121  if ( edesc_ == DESCENT_NONLINEARCG ) {
122  c2_ = p4;
123  c3_ = std::min(one-c2_,c3_);
124  }
125  }
126 
127 
128  virtual void initialize( const Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g,
130  grad_ = g.clone();
131  xtst_ = x.clone();
132  d_ = s.clone();
133  g_ = g.clone();
134  }
135 
136  virtual void run( Real &alpha, Real &fval, int &ls_neval, int &ls_ngrad,
137  const Real &gs, const Vector<Real> &s, const Vector<Real> &x,
138  Objective<Real> &obj, BoundConstraint<Real> &con ) = 0;
139 
140  // Set epsilon for epsilon active sets
141  void setData(Real &eps, const Vector<Real> &g) {
142  eps_ = eps;
143  grad_->set(g);
144  }
145 
146  // use this function to modify alpha and fval if the maximum number of iterations
147  // are reached
148  void setMaxitUpdate(Real &alpha, Real &fnew, const Real &fold) {
149  // Use local minimizer
150  if( itcond_ && acceptMin_ ) {
151  alpha = alphaMin_;
152  fnew = fmin_;
153  }
154  // Take no step
155  else if(itcond_ && !acceptMin_) {
156  alpha = 0;
157  fnew = fold;
158  }
159  setNextInitialAlpha(alpha);
160  }
161 
162 
163 protected:
164  virtual bool status( const ELineSearch type, int &ls_neval, int &ls_ngrad, const Real alpha,
165  const Real fold, const Real sgold, const Real fnew,
166  const Vector<Real> &x, const Vector<Real> &s,
167  Objective<Real> &obj, BoundConstraint<Real> &con ) {
168  Real tol = std::sqrt(ROL_EPSILON<Real>()), one(1), two(2);
169 
170  // Check Armijo Condition
171  bool armijo = false;
172  if ( con.isActivated() ) {
173  Real gs(0);
174  if ( edesc_ == DESCENT_STEEPEST ) {
175  updateIterate(*d_,x,s,alpha,con);
176  d_->scale(-one);
177  d_->plus(x);
178  gs = -s.dot(*d_);
179  }
180  else {
181  d_->set(s);
182  d_->scale(-one);
183  con.pruneActive(*d_,grad_->dual(),x,eps_);
184  gs = alpha*(grad_)->dot(d_->dual());
185  d_->zero();
186  updateIterate(*d_,x,s,alpha,con);
187  d_->scale(-one);
188  d_->plus(x);
189  con.pruneInactive(*d_,grad_->dual(),x,eps_);
190  gs += d_->dot(grad_->dual());
191  }
192  if ( fnew <= fold - c1_*gs ) {
193  armijo = true;
194  }
195  }
196  else {
197  if ( fnew <= fold + c1_*alpha*sgold ) {
198  armijo = true;
199  }
200  }
201 
202  // Check Maximum Iteration
203  itcond_ = false;
204  if ( ls_neval >= maxit_ ) {
205  itcond_ = true;
206  }
207 
208  // Check Curvature Condition
209  bool curvcond = false;
210  if ( armijo && ((type != LINESEARCH_BACKTRACKING && type != LINESEARCH_CUBICINTERP) ||
211  (edesc_ == DESCENT_NONLINEARCG)) ) {
213  if (fnew >= fold + (one-c1_)*alpha*sgold) {
214  curvcond = true;
215  }
216  }
217  else if (econd_ == CURVATURECONDITION_NULL) {
218  curvcond = true;
219  }
220  else {
221  updateIterate(*xtst_,x,s,alpha,con);
222  obj.update(*xtst_);
223  obj.gradient(*g_,*xtst_,tol);
224  Real sgnew(0);
225  if ( con.isActivated() ) {
226  d_->set(s);
227  d_->scale(-alpha);
228  con.pruneActive(*d_,s,x);
229  sgnew = -d_->dot(g_->dual());
230  }
231  else {
232  sgnew = s.dot(g_->dual());
233  }
234  ls_ngrad++;
235 
236  if ( ((econd_ == CURVATURECONDITION_WOLFE)
237  && (sgnew >= c2_*sgold))
239  && (std::abs(sgnew) <= c2_*std::abs(sgold)))
241  && (c2_*sgold <= sgnew && sgnew <= -c3_*sgold))
243  && (c2_*sgold <= sgnew && sgnew <= (two*c1_ - one)*sgold)) ) {
244  curvcond = true;
245  }
246  }
247  }
248 
249  if(fnew<fmin_) {
250  fmin_ = fnew;
251  alphaMin_ = alpha;
252  }
253 
254  if (type == LINESEARCH_BACKTRACKING || type == LINESEARCH_CUBICINTERP) {
255  if (edesc_ == DESCENT_NONLINEARCG) {
256  return ((armijo && curvcond) || itcond_);
257  }
258  else {
259  return (armijo || itcond_);
260  }
261  }
262  else {
263  return ((armijo && curvcond) || itcond_);
264  }
265  }
266 
267  virtual Real getInitialAlpha(int &ls_neval, int &ls_ngrad, const Real fval, const Real gs,
268  const Vector<Real> &x, const Vector<Real> &s,
270  Real val(1), one(1), half(0.5);
271  if (useralpha_ || usePrevAlpha_ ) {
272  val = alpha0_;
273  }
274  else {
276  Real tol = std::sqrt(ROL_EPSILON<Real>());
277  // Evaluate objective at x + s
278  updateIterate(*d_,x,s,one,con);
279  obj.update(*d_);
280  Real fnew = obj.value(*d_,tol);
281  ls_neval++;
282  // Minimize quadratic interpolate to compute new alpha
283  Real denom = (fnew - fval - gs);
284  Real alpha = ((denom > ROL_EPSILON<Real>()) ? -half*gs/denom : one);
285  val = ((alpha > alpha0bnd_) ? alpha : one);
286  }
287  else {
288  val = one;
289  }
290  }
291  return val;
292  }
293 
294  void setNextInitialAlpha( Real alpha ) {
295  if( usePrevAlpha_ ) {
296  alpha0_ = alpha;
297  }
298  }
299 
300  void updateIterate(Vector<Real> &xnew, const Vector<Real> &x, const Vector<Real> &s, Real alpha,
301  BoundConstraint<Real> &con ) {
302 
303  xnew.set(x);
304  xnew.axpy(alpha,s);
305 
306  if ( con.isActivated() ) {
307  con.project(xnew);
308  }
309 
310  }
311 
313  return itcond_ && acceptMin_;
314  }
315 
316  bool takeNoStep() {
317  return itcond_ && !acceptMin_;
318  }
319 };
320 
321 }
322 
323 #include "ROL_LineSearchFactory.hpp"
324 
325 #endif
Provides the interface to evaluate objective functions.
void setData(Real &eps, const Vector< Real > &g)
void updateIterate(Vector< Real > &xnew, const Vector< Real > &x, const Vector< Real > &s, Real alpha, BoundConstraint< Real > &con)
virtual Real getInitialAlpha(int &ls_neval, int &ls_ngrad, const Real fval, const Real gs, const Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &con)
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
bool isActivated(void) const
Check if bounds are on.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
Contains definitions of custom data types in ROL.
void pruneInactive(Vector< Real > &v, const Vector< Real > &x, Real eps=0)
Set variables to zero if they correspond to the -inactive set.
ELineSearch
Enumeration of line-search types.
Definition: ROL_Types.hpp:727
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=0)
Set variables to zero if they correspond to the -active set.
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()
ROL::Ptr< Vector< Real > > g_
EDescent StringToEDescent(std::string s)
Definition: ROL_Types.hpp:466
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Provides interface for and implements line searches.
void setMaxitUpdate(Real &alpha, Real &fnew, const Real &fold)
ROL::Ptr< Vector< Real > > grad_
Provides the interface to apply upper and lower bound constraints.
virtual ~LineSearch()
ECurvatureCondition
Enumeration of line-search curvature conditions.
Definition: ROL_Types.hpp:810
LineSearch(ROL::ParameterList &parlist)
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
ECurvatureCondition StringToECurvatureCondition(std::string s)
Definition: ROL_Types.hpp:870
virtual void run(Real &alpha, Real &fval, int &ls_neval, int &ls_ngrad, const Real &gs, const Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con)=0
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
virtual bool status(const ELineSearch type, int &ls_neval, int &ls_ngrad, const Real alpha, const Real fold, const Real sgold, const Real fnew, const Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &con)
ROL::Ptr< Vector< Real > > xtst_
ECurvatureCondition econd_
EDescent
Enumeration of descent direction types.
Definition: ROL_Types.hpp:409
ROL::Ptr< Vector< Real > > d_
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
virtual void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con)
void setNextInitialAlpha(Real alpha)