ROL
ROL_TrustRegion.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_TRUSTREGION_H
45 #define ROL_TRUSTREGION_H
46 
51 #include "ROL_Types.hpp"
52 #include "ROL_TrustRegionTypes.hpp"
53 #include "ROL_TrustRegionModel.hpp"
54 #include "ROL_ColemanLiModel.hpp"
55 #include "ROL_KelleySachsModel.hpp"
56 
57 namespace ROL {
58 
59 template<class Real>
60 class TrustRegion {
61 private:
62 
63  ROL::Ptr<Vector<Real> > prim_, dual_;
64 
66 
67  Real eta0_, eta1_, eta2_;
69  Real pRed_;
70  Real TRsafe_, eps_;
71  Real mu0_;
72 
73  std::vector<bool> useInexact_;
74 
75  Real ftol_old_;
76 
79 
80  unsigned verbosity_;
81 
82 public:
83 
84  virtual ~TrustRegion() {}
85 
86  // Constructor
87  TrustRegion( ROL::ParameterList &parlist )
88  : pRed_(0), ftol_old_(ROL_OVERFLOW<Real>()), cnt_(0), verbosity_(0) {
89  // Trust-Region Parameters
90  ROL::ParameterList list = parlist.sublist("Step").sublist("Trust Region");
91  TRmodel_ = StringToETrustRegionModel(list.get("Subproblem Model", "Kelley-Sachs"));
92  eta0_ = list.get("Step Acceptance Threshold", static_cast<Real>(0.05));
93  eta1_ = list.get("Radius Shrinking Threshold", static_cast<Real>(0.05));
94  eta2_ = list.get("Radius Growing Threshold", static_cast<Real>(0.9));
95  gamma0_ = list.get("Radius Shrinking Rate (Negative rho)", static_cast<Real>(0.0625));
96  gamma1_ = list.get("Radius Shrinking Rate (Positive rho)", static_cast<Real>(0.25));
97  gamma2_ = list.get("Radius Growing Rate", static_cast<Real>(2.5));
98  mu0_ = list.get("Sufficient Decrease Parameter", static_cast<Real>(1.e-4));
99  TRsafe_ = list.get("Safeguard Size", static_cast<Real>(100.0));
100  eps_ = TRsafe_*ROL_EPSILON<Real>();
101  // General Inexactness Information
102  ROL::ParameterList &glist = parlist.sublist("General");
103  useInexact_.clear();
104  useInexact_.push_back(glist.get("Inexact Objective Function", false));
105  useInexact_.push_back(glist.get("Inexact Gradient", false));
106  useInexact_.push_back(glist.get("Inexact Hessian-Times-A-Vector", false));
107  // Inexact Function Evaluation Information
108  ROL::ParameterList &ilist = list.sublist("Inexact").sublist("Value");
109  scale_ = ilist.get("Tolerance Scaling", static_cast<Real>(1.e-1));
110  omega_ = ilist.get("Exponent", static_cast<Real>(0.9));
111  force_ = ilist.get("Forcing Sequence Initial Value", static_cast<Real>(1.0));
112  updateIter_ = ilist.get("Forcing Sequence Update Frequency", static_cast<int>(10));
113  forceFactor_ = ilist.get("Forcing Sequence Reduction Factor", static_cast<Real>(0.1));
114  // Get verbosity level
115  verbosity_ = glist.get("Print Verbosity", 0);
116  }
117 
118  virtual void initialize( const Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g) {
119  prim_ = x.clone();
120  dual_ = g.clone();
121  }
122 
123  virtual void update( Vector<Real> &x,
124  Real &fnew,
125  Real &del,
126  int &nfval,
127  int &ngrad,
128  ETrustRegionFlag &flagTR,
129  const Vector<Real> &s,
130  const Real snorm,
131  const Real fold,
132  const Vector<Real> &g,
133  int iter,
134  Objective<Real> &obj,
136  TrustRegionModel<Real> &model ) {
137  Real tol = std::sqrt(ROL_EPSILON<Real>());
138  const Real one(1), zero(0);
139 
140  /***************************************************************************************************/
141  // BEGIN INEXACT OBJECTIVE FUNCTION COMPUTATION
142  /***************************************************************************************************/
143  // Update inexact objective function
144  Real fold1 = fold, ftol = tol; // TOL(1.e-2);
145  if ( useInexact_[0] ) {
146  if ( !(cnt_%updateIter_) && (cnt_ != 0) ) {
147  force_ *= forceFactor_;
148  }
149  //const Real oe4(1e4);
150  //Real c = scale_*std::max(TOL,std::min(one,oe4*std::max(pRed_,std::sqrt(ROL_EPSILON<Real>()))));
151  //ftol = c*std::pow(std::min(eta1_,one-eta2_)
152  // *std::min(std::max(pRed_,std::sqrt(ROL_EPSILON<Real>())),force_),one/omega_);
153  //if ( ftol_old_ > ftol || cnt_ == 0 ) {
154  // ftol_old_ = ftol;
155  // fold1 = obj.value(x,ftol_old_);
156  //}
157  //cnt_++;
158  Real eta = static_cast<Real>(0.999)*std::min(eta1_,one-eta2_);
159  ftol = scale_*std::pow(eta*std::min(pRed_,force_),one/omega_);
160  ftol_old_ = ftol;
161  fold1 = obj.value(x,ftol_old_);
162  cnt_++;
163  }
164  // Evaluate objective function at new iterate
165  prim_->set(x); prim_->plus(s);
166  obj.update(*prim_,true);
167  fnew = obj.value(*prim_,ftol);
168 
169  nfval = 1;
170  Real aRed = fold1 - fnew;
171  /***************************************************************************************************/
172  // FINISH INEXACT OBJECTIVE FUNCTION COMPUTATION
173  /***************************************************************************************************/
174 
175  /***************************************************************************************************/
176  // BEGIN COMPUTE RATIO OF ACTUAL AND PREDICTED REDUCTION
177  /***************************************************************************************************/
178  // Modify Actual and Predicted Reduction According to Model
179  model.updateActualReduction(aRed,s);
181 
182  if ( verbosity_ > 0 ) {
183  std::cout << std::endl;
184  std::cout << " Computation of actual and predicted reduction" << std::endl;
185  std::cout << " Current objective function value: " << fold1 << std::endl;
186  std::cout << " New objective function value: " << fnew << std::endl;
187  std::cout << " Actual reduction: " << aRed << std::endl;
188  std::cout << " Predicted reduction: " << pRed_ << std::endl;
189  }
190 
191  // Compute Ratio of Actual and Predicted Reduction
192  Real EPS = eps_*((one > std::abs(fold1)) ? one : std::abs(fold1));
193  Real aRed_safe = aRed + EPS, pRed_safe = pRed_ + EPS;
194  Real rho(0);
195  if (((std::abs(aRed_safe) < eps_) && (std::abs(pRed_safe) < eps_)) || aRed == pRed_) {
196  rho = one;
197  flagTR = TRUSTREGION_FLAG_SUCCESS;
198  }
199  else if ( std::isnan(aRed_safe) || std::isnan(pRed_safe) ) {
200  rho = -one;
201  flagTR = TRUSTREGION_FLAG_NAN;
202  }
203  else {
204  rho = aRed_safe/pRed_safe;
205  if (pRed_safe < zero && aRed_safe > zero) {
207  }
208  else if (aRed_safe <= zero && pRed_safe > zero) {
210  }
211  else if (aRed_safe <= zero && pRed_safe < zero) {
213  }
214  else {
215  flagTR = TRUSTREGION_FLAG_SUCCESS;
216  }
217  }
218 
219  if ( verbosity_ > 0 ) {
220  std::cout << " Safeguard: " << eps_ << std::endl;
221  std::cout << " Actual reduction with safeguard: " << aRed_safe << std::endl;
222  std::cout << " Predicted reduction with safeguard: " << pRed_safe << std::endl;
223  std::cout << " Ratio of actual and predicted reduction: " << rho << std::endl;
224  std::cout << " Trust-region flag: " << flagTR << std::endl;
225  }
226  /***************************************************************************************************/
227  // FINISH COMPUTE RATIO OF ACTUAL AND PREDICTED REDUCTION
228  /***************************************************************************************************/
229 
230 
231  /***************************************************************************************************/
232  // BEGIN CHECK SUFFICIENT DECREASE FOR BOUND CONSTRAINED PROBLEMS
233  /***************************************************************************************************/
234  bool decr = true;
236  if ( rho >= eta0_ && (std::abs(aRed_safe) > eps_) ) {
237  // Compute Criticality Measure || x - P( x - g ) ||
238  prim_->set(x);
239  prim_->axpy(-one,g.dual());
240  bnd.project(*prim_);
241  prim_->scale(-one);
242  prim_->plus(x);
243  Real pgnorm = prim_->norm();
244  // Compute Scaled Measure || x - P( x - lam * PI(g) ) ||
245  prim_->set(g.dual());
246  bnd.pruneActive(*prim_,g,x);
247  Real lam = std::min(one, del/prim_->norm());
248  prim_->scale(-lam);
249  prim_->plus(x);
250  bnd.project(*prim_);
251  prim_->scale(-one);
252  prim_->plus(x);
253  pgnorm *= prim_->norm();
254  // Sufficient decrease?
255  decr = ( aRed_safe >= mu0_*pgnorm );
256  flagTR = (!decr ? TRUSTREGION_FLAG_QMINSUFDEC : flagTR);
257 
258  if ( verbosity_ > 0 ) {
259  std::cout << " Decrease lower bound (constraints): " << mu0_*pgnorm << std::endl;
260  std::cout << " Trust-region flag (constraints): " << flagTR << std::endl;
261  std::cout << " Is step feasible: " << bnd.isFeasible(x) << std::endl;
262  }
263  }
264  }
265  /***************************************************************************************************/
266  // FINISH CHECK SUFFICIENT DECREASE FOR BOUND CONSTRAINED PROBLEMS
267  /***************************************************************************************************/
268 
269  /***************************************************************************************************/
270  // BEGIN STEP ACCEPTANCE AND TRUST REGION RADIUS UPDATE
271  /***************************************************************************************************/
272  if ( verbosity_ > 0 ) {
273  std::cout << " Norm of step: " << snorm << std::endl;
274  std::cout << " Trust-region radius before update: " << del << std::endl;
275  }
276  if ((rho < eta0_ && flagTR == TRUSTREGION_FLAG_SUCCESS) || flagTR >= 2 || !decr ) { // Step Rejected
277  fnew = fold1;
278  if (rho < zero) { // Negative reduction, interpolate to find new trust-region radius
279  Real gs(0);
280  if ( bnd.isActivated() ) {
281  model.dualTransform(*dual_, *model.getGradient());
282  gs = dual_->dot(s.dual());
283  }
284  else {
285  gs = g.dot(s.dual());
286  }
287  Real modelVal = model.value(s,tol);
288  modelVal += fold1;
289  Real theta = (one-eta2_)*gs/((one-eta2_)*(fold1+gs)+eta2_*modelVal-fnew);
290  del = std::min(gamma1_*std::min(snorm,del),std::max(gamma0_,theta)*del);
291  if ( verbosity_ > 0 ) {
292  std::cout << " Interpolation model value: " << modelVal << std::endl;
293  std::cout << " Interpolation step length: " << theta << std::endl;
294  }
295  }
296  else { // Shrink trust-region radius
297  del = gamma1_*std::min(snorm,del);
298  }
299  obj.update(x,true,iter);
300  }
301  else if ((rho >= eta0_ && flagTR != TRUSTREGION_FLAG_NPOSPREDNEG) ||
302  (flagTR == TRUSTREGION_FLAG_POSPREDNEG)) { // Step Accepted
303  x.plus(s);
304  obj.update(x,true,iter);
305  if (rho >= eta2_) { // Increase trust-region radius
306  del = gamma2_*del;
307  }
308  }
309  if ( verbosity_ > 0 ) {
310  std::cout << " Trust-region radius after update: " << del << std::endl;
311  std::cout << std::endl;
312  }
313  /***************************************************************************************************/
314  // FINISH STEP ACCEPTANCE AND TRUST REGION RADIUS UPDATE
315  /***************************************************************************************************/
316  }
317 
318  virtual void run( Vector<Real> &s, // Step (to be computed)
319  Real &snorm, // Step norm (to be computed)
320  int &iflag, // Exit flag (to be computed)
321  int &iter, // Iteration count (to be computed)
322  const Real del, // Trust-region radius
323  TrustRegionModel<Real> &model ) = 0; // Trust-region model
324 
325  void setPredictedReduction(const Real pRed) {
326  pRed_ = pRed;
327  }
328 
329  Real getPredictedReduction(void) const {
330  return pRed_;
331  }
332 };
333 
334 }
335 
337 
338 #endif
Provides the interface to evaluate objective functions.
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 ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void plus(const Vector &x)=0
Compute , where .
bool isActivated(void) const
Check if bounds are on.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
virtual void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
Contains definitions of custom data types in ROL.
Provides interface for and implements trust-region subproblem solvers.
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=0)
Set variables to zero if they correspond to the -active set.
Provides the interface to evaluate trust-region model functions.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
virtual const Ptr< const Vector< Real > > getGradient(void) const
virtual void updatePredictedReduction(Real &pred, const Vector< Real > &s)
virtual Real dot(const Vector &x) const =0
Compute where .
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
Real getPredictedReduction(void) const
ETrustRegionModel StringToETrustRegionModel(std::string s)
virtual Real value(const Vector< Real > &s, Real &tol)
Compute value.
void setPredictedReduction(const Real pRed)
Real ROL_OVERFLOW(void)
Platform-dependent maximum double.
Definition: ROL_Types.hpp:102
ETrustRegionModel
Enumeration of trust-region model types.
TrustRegion(ROL::ParameterList &parlist)
Provides the interface to apply upper and lower bound constraints.
virtual void updateActualReduction(Real &ared, const Vector< Real > &s)
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
ETrustRegionFlag
Enumation of flags used by trust-region solvers.
ETrustRegionModel TRmodel_
ROL::Ptr< Vector< Real > > dual_
virtual bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
Contains definitions of enums for trust region algorithms.
std::vector< bool > useInexact_
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
virtual void run(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)=0
virtual void update(Vector< Real > &x, Real &fnew, Real &del, int &nfval, int &ngrad, ETrustRegionFlag &flagTR, const Vector< Real > &s, const Real snorm, const Real fold, const Vector< Real > &g, int iter, Objective< Real > &obj, BoundConstraint< Real > &bnd, TrustRegionModel< Real > &model)
ROL::Ptr< Vector< Real > > prim_