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 "Teuchos_RCP.hpp"
52 #include "Teuchos_ParameterList.hpp"
53 #include "ROL_Types.hpp"
54 #include "ROL_Vector.hpp"
55 #include "ROL_Objective.hpp"
56 #include "ROL_BoundConstraint.hpp"
57 
58 namespace ROL {
59 
60 template<class Real>
61 class LineSearch {
62 private:
63 
66 
67  bool useralpha_;
68  Real alpha0_;
69  int maxit_;
70  Real c1_;
71  Real c2_;
72  Real c3_;
73  Real eps_;
74 
75  Teuchos::RCP<Vector<Real> > xtst_;
76  Teuchos::RCP<Vector<Real> > d_;
77  Teuchos::RCP<Vector<Real> > g_;
78  Teuchos::RCP<const Vector<Real> > grad_;
79 
80 public:
81 
82  virtual ~LineSearch() {}
83 
84  // Constructor
85  LineSearch( Teuchos::ParameterList &parlist ) : eps_(0.0) {
86  // Enumerations
87  edesc_ = StringToEDescent(parlist.get("Descent Type","Quasi-Newton Method"));
88  econd_ = StringToECurvatureCondition( parlist.get("Linesearch Curvature Condition","Strong Wolfe Conditions"));
89  // Linesearc Parameters
90  maxit_ = parlist.get("Maximum Number of Function Evaluations", 20);
91  c1_ = parlist.get("Sufficient Decrease Parameter", 1.e-4);
92  c2_ = parlist.get("Curvature Conditions Parameter", 0.9);
93  c3_ = parlist.get("Curvature Conditions Parameter: Generalized Wolfe", 0.6);
94  alpha0_ = parlist.get("Initial Linesearch Parameter",1.0);
95  useralpha_ = parlist.get("User Defined Linesearch Parameter",false);
96 
97  if ( c1_ < 0.0 ) {
98  c1_ = 1.e-4;
99  }
100  if ( c2_ < 0.0 ) {
101  c2_ = 0.9;
102  }
103  if ( c3_ < 0.0 ) {
104  c3_ = 0.9;
105  }
106  if ( c2_ <= c1_ ) {
107  c1_ = 1.e-4;
108  c2_ = 0.9;
109  }
110  if ( edesc_ == DESCENT_NONLINEARCG ) {
111  c2_ = 0.4;
112  c3_ = std::min(1.0-c2_,c3_);
113  }
114  }
115 
116  virtual void initialize( const Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g,
118  grad_ = Teuchos::rcp(&g, false);
119  xtst_ = x.clone();
120  d_ = s.clone();
121  g_ = g.clone();
122  }
123 
124  void setData(Real &eps) {
125  eps_ = eps;
126  }
127 
128  virtual bool status( const ELineSearch type, int &ls_neval, int &ls_ngrad, const Real alpha,
129  const Real fold, const Real sgold, const Real fnew,
130  const Vector<Real> &x, const Vector<Real> &s,
131  Objective<Real> &obj, BoundConstraint<Real> &con ) {
132  Real tol = std::sqrt(ROL_EPSILON);
133 
134  // Check Armijo Condition
135  bool armijo = false;
136  if ( con.isActivated() ) {
137  Real gs = 0.0;
138  if ( edesc_ == DESCENT_STEEPEST ) {
139  updateIterate(*d_,x,s,alpha,con);
140  d_->scale(-1.0);
141  d_->plus(x);
142  gs = -s.dot(*d_);
143  }
144  else {
145  d_->set(s);
146  d_->scale(-1.0);
147  con.pruneActive(*d_,*(grad_),x,eps_);
148  gs = alpha*(grad_)->dot(*d_);
149  d_->zero();
150  updateIterate(*d_,x,s,alpha,con);
151  d_->scale(-1.0);
152  d_->plus(x);
153  con.pruneInactive(*d_,*(grad_),x,eps_);
154  gs += d_->dot(grad_->dual());
155  }
156  if ( fnew <= fold - c1_*gs ) {
157  armijo = true;
158  }
159  }
160  else {
161  if ( fnew <= fold + c1_*alpha*sgold ) {
162  armijo = true;
163  }
164  }
165 
166  // Check Maximum Iteration
167  bool itcond = false;
168  if ( ls_neval >= maxit_ ) {
169  itcond = true;
170  }
171 
172  // Check Curvature Condition
173  bool curvcond = false;
174  if ( armijo && ((type != LINESEARCH_BACKTRACKING && type != LINESEARCH_CUBICINTERP) ||
175  (edesc_ == DESCENT_NONLINEARCG)) ) {
177  if (fnew >= fold + (1.0-c1_)*alpha*sgold) {
178  curvcond = true;
179  }
180  }
181  else if (econd_ == CURVATURECONDITION_NULL) {
182  curvcond = true;
183  }
184  else {
185  updateIterate(*xtst_,x,s,alpha,con);
186  obj.update(*xtst_);
187  obj.gradient(*g_,*xtst_,tol);
188  Real sgnew = 0.0;
189  if ( con.isActivated() ) {
190  d_->set(s);
191  d_->scale(-alpha);
192  con.pruneActive(*d_,s,x);
193  sgnew = -d_->dot(g_->dual());
194  }
195  else {
196  sgnew = s.dot(g_->dual());
197  }
198  ls_ngrad++;
199 
200  if ( ((econd_ == CURVATURECONDITION_WOLFE)
201  && (sgnew >= c2_*sgold))
203  && (std::abs(sgnew) <= c2_*std::abs(sgold)))
205  && (c2_*sgold <= sgnew && sgnew <= -c3_*sgold))
207  && (c2_*sgold <= sgnew && sgnew <= (2.0*c1_ - 1.0)*sgold)) ) {
208  curvcond = true;
209  }
210  }
211  }
212 
213  if (type == LINESEARCH_BACKTRACKING || type == LINESEARCH_CUBICINTERP) {
214  if (edesc_ == DESCENT_NONLINEARCG) {
215  return ((armijo && curvcond) || itcond);
216  }
217  else {
218  return (armijo || itcond);
219  }
220  }
221  else {
222  return ((armijo && curvcond) || itcond);
223  }
224  }
225 
226  virtual void run( Real &alpha, Real &fval, int &ls_neval, int &ls_ngrad,
227  const Real &gs, const Vector<Real> &s, const Vector<Real> &x,
228  Objective<Real> &obj, BoundConstraint<Real> &con ) = 0;
229 
230  virtual Real getInitialAlpha(int &ls_neval, int &ls_ngrad, const Real fval, const Real gs,
231  const Vector<Real> &x, const Vector<Real> &s,
233  Real val = 1.0;
234  if (useralpha_) {
235  val = alpha0_;
236  }
237  else {
239  Real tol = std::sqrt(ROL_EPSILON);
240  Real alpha = 1.0;
241  // Evaluate objective at x + s
242  updateIterate(*d_,x,s,alpha,con);
243  obj.update(*d_);
244  Real fnew = obj.value(*d_,tol);
245  ls_neval++;
246  // Minimize quadratic interpolate to compute new alpha
247  alpha = -gs/(2.0*(fnew-fval-gs));
248  // Evaluate objective at x + alpha s
249  updateIterate(*d_,x,s,alpha,con);
250  obj.update(*d_);
251  fnew = obj.value(*d_,tol);
252  ls_neval++;
253  // Ensure that sufficient decrease and curvature conditions are satisfied
254  bool stat = status(LINESEARCH_BISECTION,ls_neval,ls_ngrad,alpha,fval,gs,fnew,x,s,obj,con);
255  if ( !stat ) {
256  alpha = 1.0;
257  }
258  val = alpha;
259  }
260  else {
261  val = 1.0;
262  }
263  }
264  return val;
265  }
266 
267  void updateIterate(Vector<Real> &xnew, const Vector<Real> &x, const Vector<Real> &s, Real alpha,
268  BoundConstraint<Real> &con ) {
269  xnew.set(x);
270  xnew.axpy(alpha,s);
271  if ( con.isActivated() ) {
272  con.project(xnew);
273  }
274  }
275 };
276 
277 }
278 
279 #include "ROL_IterationScaling.hpp"
281 #include "ROL_BackTracking.hpp"
282 #include "ROL_CubicInterp.hpp"
283 #include "ROL_Bisection.hpp"
284 #include "ROL_GoldenSection.hpp"
285 #include "ROL_Brents.hpp"
286 
287 #endif
Provides the interface to evaluate objective functions.
void updateIterate(Vector< Real > &xnew, const Vector< Real > &x, const Vector< Real > &s, Real alpha, BoundConstraint< Real > &con)
void setData(Real &eps)
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 void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:141
Teuchos::RCP< Vector< Real > > xtst_
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
Contains definitions of custom data types in ROL.
Teuchos::RCP< Vector< Real > > g_
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
ELineSearch
Enumeration of line-search types.
Definition: ROL_Types.hpp:527
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
LineSearch(Teuchos::ParameterList &parlist)
virtual Real dot(const Vector &x) const =0
Compute where .
EDescent StringToEDescent(std::string s)
Definition: ROL_Types.hpp:277
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
bool isActivated(void)
Check if bounds are on.
Provides interface for and implements line searches.
void pruneInactive(Vector< Real > &v, const Vector< Real > &x, Real eps=0.0)
Set variables to zero if they correspond to the -inactive set.
Teuchos::RCP< Vector< Real > > d_
Provides the interface to apply upper and lower bound constraints.
virtual ~LineSearch()
ECurvatureCondition
Enumeration of line-search curvature conditions.
Definition: ROL_Types.hpp:610
Teuchos::RCP< const Vector< Real > > grad_
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:194
virtual void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=0.0)
Set variables to zero if they correspond to the -active set.
ECurvatureCondition StringToECurvatureCondition(std::string s)
Definition: ROL_Types.hpp:670
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)
ECurvatureCondition econd_
EDescent
Enumeration of descent direction types.
Definition: ROL_Types.hpp:220
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)
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115