ROL
ROL_BundleStep.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_BUNDLE_STEP_H
45 #define ROL_BUNDLE_STEP_H
46 
47 #include "ROL_Bundle.hpp"
48 #include "ROL_Types.hpp"
49 #include "ROL_Step.hpp"
50 #include "ROL_Vector.hpp"
51 #include "ROL_Objective.hpp"
52 #include "ROL_BoundConstraint.hpp"
53 #include "ROL_LineSearch.hpp"
54 
55 #include "Teuchos_ParameterList.hpp"
56 #include "Teuchos_RCP.hpp"
57 
63 namespace ROL {
64 
65 template <class Real>
66 class BundleStep : public Step<Real> {
67 private:
68  // Bundle
69  Teuchos::RCP<Bundle<Real> > bundle_; // Bundle of subgradients and linearization errors
70  Teuchos::RCP<LineSearch<Real> > lineSearch_; // Line-search object for nonconvex problems
71 
72  // Dual cutting plane solution
73  unsigned QPiter_; // Number of QP solver iterations
74  unsigned QPmaxit_; // Maximum number of QP iterations
75  Real QPtol_; // QP subproblem tolerance
76 
77  // Step flag
78  int step_flag_; // Whether serious or null step
79 
80  // Additional storage
81  Teuchos::RCP<Vector<Real> > y_;
82 
83  // Updated iterate storage
84  Real linErrNew_;
85  Real valueNew_;
86 
87  // Aggregate subgradients, linearizations, and distance measures
88  Teuchos::RCP<Vector<Real> > aggSubGradNew_; // New aggregate subgradient
89  Real aggSubGradOldNorm_; // Old aggregate subgradient norm
90  Real aggLinErrNew_; // New aggregate linearization error
91  Real aggLinErrOld_; // Old aggregate linearization error
92  Real aggDistMeasNew_; // New aggregate distance measure
93 
94  // Algorithmic parameters
95  Real T_;
96  Real tol_;
97  Real m1_;
98  Real m2_;
99  Real m3_;
100  Real nu_;
101 
102  // Line-search parameters
104 
106  bool isConvex_;
107 
108  Real ftol_;
109 
110  void LineSearchFactory(Teuchos::ParameterList &parlist) {
111  ls_maxit_ = parlist.get("Maximum Number of Function Evaluations",20);
112  ELineSearch els = StringToELineSearch(parlist.get("Linesearch Type","Cubic Interpolation"));
113  switch(els) {
115  lineSearch_ = Teuchos::rcp( new IterationScaling<Real>(parlist) ); break;
117  lineSearch_ = Teuchos::rcp( new PathBasedTargetLevel<Real>(parlist) ); break;
119  lineSearch_ = Teuchos::rcp( new BackTracking<Real>(parlist) ); break;
121  lineSearch_ = Teuchos::rcp( new Bisection<Real>(parlist) ); break;
122  case LINESEARCH_BRENTS:
123  lineSearch_ = Teuchos::rcp( new Brents<Real>(parlist) ); break;
125  lineSearch_ = Teuchos::rcp( new GoldenSection<Real>(parlist) ); break;
127  default:
128  lineSearch_ = Teuchos::rcp( new CubicInterp<Real>(parlist) ); break;
129  }
130  }
131 
132 public:
133 
134  BundleStep(Teuchos::ParameterList &parlist)
135  : QPiter_(0), step_flag_(0), aggLinErrNew_(0.0), aggLinErrOld_(0.0), aggDistMeasNew_(0.0),
136  first_print_(true), ftol_(ROL_EPSILON) {
137  Teuchos::RCP<StepState<Real> > state = Step<Real>::getState();
138  state->searchSize = parlist.get("Bundle Step: Initial Trust-Region Parameter",1.e3);
139  T_ = parlist.get("Bundle Step: Maximum Trust-Region Parameter",1.e8);
140  tol_ = parlist.get("Bundle Step: Epsilon Solution Tolerance",1.e-6);
141  m1_ = parlist.get("Bundle Step: Upper Threshold for Serious Step",0.1);
142  m2_ = parlist.get("Bundle Step: Lower Threshold for Serious Step",0.2);
143  m3_ = parlist.get("Bundle Step: Upper Threshold for Null Step",0.9);
144  nu_ = parlist.get("Bundle Step: Tolerance for Trust-Region Parameter",1.e-3);
145 
146  // Initialize bundle
147  Real coeff = parlist.get("Bundle Step: Distance Measure Coefficient",0.0);
148  unsigned maxSize = parlist.get("Bundle Step: Maximum Bundle Size",200);
149  unsigned remSize = parlist.get("Bundle Step: Removal Size for Bundle Update",2);
150  bundle_ = Teuchos::rcp(new Bundle<Real>(maxSize,coeff,remSize));
151  isConvex_ = ((coeff == 0.0) ? true : false);
152 
153  // Initialize QP solver
154  QPtol_ = parlist.get("Bundle Step: Cutting Plane Solver Tolerance",1.e-8);
155  QPmaxit_ = parlist.get("Bundle Step: Cutting Plane Solver Maximum Number of Iterations",1000);
156 
157  // Initialize Line Search
158  lineSearch_ = Teuchos::null;
159  if ( !isConvex_ ) {
160  LineSearchFactory(parlist);
161  }
162  }
163 
164  void initialize( Vector<Real> &x, const Vector<Real> &g,
166  AlgorithmState<Real> &algo_state ) {
167  // Call default initializer, but maintain current searchSize
168  Teuchos::RCP<StepState<Real> > state = Step<Real>::getState();
169  Real searchSize = state->searchSize;
170  Step<Real>::initialize(x,x,g,obj,con,algo_state);
171  state->searchSize = searchSize;
172  // Initialize bundle
173  bundle_->initialize(*(state->gradientVec));
174  // Initialize storage for updated iterate
175  y_ = x.clone();
176  // Initialize storage for aggregate subgradients
177  aggSubGradNew_ = x.clone();
178  aggSubGradOldNorm_ = algo_state.gnorm;
179  }
180 
182  BoundConstraint<Real> &con, AlgorithmState<Real> &algo_state ) {
183  Teuchos::RCP<StepState<Real> > state = Step<Real>::getState();
184  first_print_ = false; // Print header only on first serious step
185  QPiter_ = (step_flag_ ? 0 : QPiter_); // Reset QPiter only on serious steps
186  Real v = 0.0, l = 0.0, u = T_, gd = 0.0; // Scalar storage
187  bool flag = true;
188  while (flag) {
189  /*************************************************************/
190  /******** Solve Dual Cutting Plane QP Problem ****************/
191  /*************************************************************/
192  QPiter_ += bundle_->solveDual(state->searchSize,QPmaxit_,QPtol_); // Solve QP subproblem
193  bundle_->aggregate(*aggSubGradNew_,aggLinErrNew_,aggDistMeasNew_); // Compute aggregate info
194  algo_state.aggregateGradientNorm = aggSubGradNew_->norm(); // Aggregate subgradient norm
195  /*************************************************************/
196  /******** Construct Cutting Plane Solution *******************/
197  /*************************************************************/
198  v = -state->searchSize*std::pow(algo_state.aggregateGradientNorm,2.0)-aggLinErrNew_; // CP objective value
199  s.set(*aggSubGradNew_); s.scale(-state->searchSize); // CP solution
200  algo_state.snorm = state->searchSize*algo_state.aggregateGradientNorm; // Step norm
201  /*************************************************************/
202  /******** Decide Whether Step is Serious or Null *************/
203  /*************************************************************/
204  if (std::max(algo_state.aggregateGradientNorm,aggLinErrNew_) <= tol_) {
205  // Current iterate is already epsilon optimal!
206  s.zero(); algo_state.snorm = 0.0;
207  flag = false;
208  step_flag_ = 1;
209  algo_state.flag = true;
210  break;
211  }
212  else {
213  // Current iterate is not epsilon optimal.
214  y_->set(x); y_->plus(s); // y is the candidate iterate
215  obj.update(*y_,true,algo_state.iter); // Update objective at y
216  valueNew_ = obj.value(*y_,ftol_); // Compute objective value at y
217  algo_state.nfval++;
218  obj.gradient(*(state->gradientVec),*y_,ftol_); // Compute objective (sub)gradient at y
219  algo_state.ngrad++;
220  // Compute new linearization error and distance measure
221  gd = s.dot(*(state->gradientVec));
222  linErrNew_ = algo_state.value - (valueNew_ - gd); // Linearization error
223  // Determine whether to take a serious or null step
224  bool SS1 = (valueNew_-algo_state.value < m1_*v);
225  //bool NS1 = (valueNew_-algo_state.value >= m1_*v);
226  bool NS2a = (bundle_->computeAlpha(algo_state.snorm,linErrNew_) <= m3_*aggLinErrOld_);
227  bool NS2b = (std::abs(algo_state.value-valueNew_) <= aggSubGradOldNorm_ + aggLinErrOld_);
228  if ( isConvex_ ) {
229  /************* Convex objective ****************/
230  if ( SS1 ) {
231  bool SS2 = (gd >= m2_*v || state->searchSize >= T_-nu_);
232  if ( SS2 ) { // Serious Step
233  step_flag_ = 1;
234  flag = false;
235  break;
236  }
237  else { // Increase trust-region radius
238  l = state->searchSize;
239  state->searchSize = 0.5*(u+l);
240  }
241  }
242  else {
243  if ( NS2a || NS2b ) { // Null step
244  s.zero(); algo_state.snorm = 0.0;
245  step_flag_ = 0;
246  flag = false;
247  break;
248  }
249  else { // Decrease trust-region radius
250  u = state->searchSize;
251  state->searchSize = 0.5*(u+l);
252  }
253  }
254  }
255  else {
256  /************* Nonconvex objective *************/
257  bool NS3 = (gd - bundle_->computeAlpha(algo_state.snorm,linErrNew_) >= m2_*v);
258  if ( SS1 ) { // Serious step
259  step_flag_ = 1;
260  flag = false;
261  break;
262  }
263  else {
264  if ( NS2a || NS2b ) {
265  if ( NS3 ) { // Null step
266  s.zero();
267  step_flag_ = 0;
268  flag = false;
269  break;
270  }
271  else {
272  if ( NS2b ) { // Line search
273  Real alpha = 0.0;
274  int ls_nfval = 0, ls_ngrad = 0;
275  lineSearch_->run(alpha,valueNew_,ls_nfval,ls_ngrad,gd,s,x,obj,con);
276  if ( ls_nfval == ls_maxit_ ) { // Null step
277  s.zero();
278  step_flag_ = 0;
279  flag = false;
280  break;
281  }
282  else { // Serious step
283  s.scale(alpha);
284  step_flag_ = 1;
285  flag = false;
286  break;
287  }
288  }
289  else { // Decrease trust-region radius
290  u = state->searchSize;
291  state->searchSize = 0.5*(u+l);
292  }
293  }
294  }
295  else { // Decrease trust-region radius
296  u = state->searchSize;
297  state->searchSize = 0.5*(u+l);
298  }
299  }
300  }
301  }
302  } // End While
303  /*************************************************************/
304  /******** Update Algorithm State *****************************/
305  /*************************************************************/
306  algo_state.aggregateModelError = aggLinErrNew_;
309  } // End Compute
310 
311  void update( Vector<Real> &x, const Vector<Real> &s, Objective<Real> &obj,
312  BoundConstraint<Real> &con, AlgorithmState<Real> &algo_state ) {
313  Teuchos::RCP<StepState<Real> > state = Step<Real>::getState();
314  if ( !algo_state.flag ) {
315  /*************************************************************/
316  /******** Reset Bundle If Maximum Size Reached ***************/
317  /*************************************************************/
318  bundle_->reset(*aggSubGradNew_,aggLinErrNew_,algo_state.snorm);
319  /*************************************************************/
320  /******** Update Bundle with Step Information ****************/
321  /*************************************************************/
322  if ( step_flag_ ) {
323  // Serious step was taken
324  x.plus(s); // Update iterate
325  Real valueOld = algo_state.value; // Store previous objective value
326  algo_state.value = valueNew_; // Add new objective value to state
327  bundle_->update(step_flag_,valueNew_-valueOld,algo_state.snorm,*(state->gradientVec),s);
328  }
329  else {
330  // Null step was taken
331  bundle_->update(step_flag_,linErrNew_,algo_state.snorm,*(state->gradientVec),s);
332  }
333  }
334  /*************************************************************/
335  /******** Update Algorithm State *****************************/
336  /*************************************************************/
337  algo_state.iterateVec->set(x);
338  algo_state.gnorm = (state->gradientVec)->norm();
339  if ( step_flag_ ) {
340  algo_state.iter++;
341  }
342  } // End Update
343 
344  std::string printHeader( void ) const {
345  std::stringstream hist;
346  hist << " ";
347  hist << std::setw(6) << std::left << "iter";
348  hist << std::setw(15) << std::left << "value";
349  hist << std::setw(15) << std::left << "gnorm";
350  hist << std::setw(15) << std::left << "snorm";
351  hist << std::setw(10) << std::left << "#fval";
352  hist << std::setw(10) << std::left << "#grad";
353  hist << std::setw(15) << std::left << "znorm";
354  hist << std::setw(15) << std::left << "alpha";
355  hist << std::setw(15) << std::left << "TRparam";
356  hist << std::setw(10) << std::left << "QPiter";
357  hist << "\n";
358  return hist.str();
359  }
360 
361  std::string printName( void ) const {
362  std::stringstream hist;
363  hist << "\n" << "Bundle Trust-Region Algorithm \n";
364  return hist.str();
365  }
366 
367  std::string print( AlgorithmState<Real> &algo_state, bool print_header = false ) const {
368  Teuchos::RCP<const StepState<Real> > state = Step<Real>::getStepState();
369  std::stringstream hist;
370  hist << std::scientific << std::setprecision(6);
371  if ( algo_state.iter == 0 && first_print_ ) {
372  hist << printName();
373  if ( print_header ) {
374  hist << printHeader();
375  }
376  hist << " ";
377  hist << std::setw(6) << std::left << algo_state.iter;
378  hist << std::setw(15) << std::left << algo_state.value;
379  hist << std::setw(15) << std::left << algo_state.gnorm;
380  hist << "\n";
381  }
382  if ( step_flag_ && algo_state.iter > 0 ) {
383  if ( print_header ) {
384  hist << printHeader();
385  }
386  else {
387  hist << " ";
388  hist << std::setw(6) << std::left << algo_state.iter;
389  hist << std::setw(15) << std::left << algo_state.value;
390  hist << std::setw(15) << std::left << algo_state.gnorm;
391  hist << std::setw(15) << std::left << algo_state.snorm;
392  hist << std::setw(10) << std::left << algo_state.nfval;
393  hist << std::setw(10) << std::left << algo_state.ngrad;
394  hist << std::setw(15) << std::left << algo_state.aggregateGradientNorm;
395  hist << std::setw(15) << std::left << algo_state.aggregateModelError;
396  hist << std::setw(15) << std::left << state->searchSize;
397  hist << std::setw(10) << std::left << QPiter_;
398  hist << "\n";
399  }
400  }
401  return hist.str();
402  };
403 
404 }; // class BundleStep
405 
406 } // namespace ROL
407 
408 #endif
Provides the interface to evaluate objective functions.
virtual void scale(const Real alpha)=0
Compute where .
virtual void plus(const Vector &x)=0
Compute , where .
ELineSearch StringToELineSearch(std::string s)
Definition: ROL_Types.hpp:593
Implements a golden section line search.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
Provides the interface to compute optimization steps.
Definition: ROL_Step.hpp:63
Teuchos::RCP< StepState< Real > > getState(void)
Definition: ROL_Step.hpp:68
Contains definitions of custom data types in ROL.
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
std::string printName(void) const
Print step name.
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:155
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
virtual Real dot(const Vector &x) const =0
Compute where .
State for algorithm class. Will be used for restarts.
Definition: ROL_Types.hpp:76
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Provides an implementation of path-based target leve line search.
virtual Teuchos::RCP< const StepState< Real > > getStepState(void) const
Get state for step object.
Definition: ROL_Step.hpp:188
Implements a Brent's method line search.
Definition: ROL_Brents.hpp:57
void LineSearchFactory(Teuchos::ParameterList &parlist)
Implements cubic interpolation back tracking line search.
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Initialize step with bound constraint.
Teuchos::RCP< Vector< Real > > y_
Implements a simple back tracking line search.
std::string printHeader(void) const
Print iterate header.
Provides the interface to apply upper and lower bound constraints.
void update(Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Update step, if successful.
std::string print(AlgorithmState< Real > &algo_state, bool print_header=false) const
Print iterate status.
void compute(Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Compute step.
Provides the interface to compute bundle trust-region steps.
Teuchos::RCP< Bundle< Real > > bundle_
virtual void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Initialize step with bound constraint.
Definition: ROL_Step.hpp:83
Teuchos::RCP< Vector< Real > > iterateVec
Definition: ROL_Types.hpp:90
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:194
Implements a bisection line search.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Provides an implementation of iteration scaled line search.
Teuchos::RCP< Vector< Real > > aggSubGradNew_
Teuchos::RCP< LineSearch< Real > > lineSearch_
BundleStep(Teuchos::ParameterList &parlist)
Provides the interface for and implments a bundle.
Definition: ROL_Bundle.hpp:62
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115