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_HelperFunctions.hpp"
53 
54 namespace ROL {
55 
56 template<class Real>
57 class TrustRegion {
58 private:
59 
60  Teuchos::RCP<Vector<Real> > xupdate_;
61  Teuchos::RCP<Vector<Real> > Hs_;
62 
63  Real delmax_;
64  Real eta0_;
65  Real eta1_;
66  Real eta2_;
67  Real gamma0_;
68  Real gamma1_;
69  Real gamma2_;
70 
71  Real pRed_;
72 
73  Real TRsafe_;
74  Real eps_;
75 
76  std::vector<bool> useInexact_;
77 
78  Real ftol_old_;
79 
80  Real scale_;
81  Real omega_;
82  Real force_;
85  int cnt_;
86 
87  bool softUp_;
88 
89  void updateObj( Vector<Real> &x, int iter, ProjectedObjective<Real> &pObj ) {
90  if ( !softUp_ ) {
91  pObj.update(x,true,iter);
92  }
93  else {
94  pObj.update(x);
95  }
96  }
97 
98 
99 public:
100 
101  virtual ~TrustRegion() {}
102 
103  // Constructor
104  TrustRegion( Teuchos::ParameterList & parlist ) : ftol_old_(ROL_OVERFLOW), cnt_(0) {
105  // Unravel Parameter List
106  // Trust-Region Parameters
107  delmax_ = parlist.get("Maximum Trust-Region Radius", 5000.0);
108  eta0_ = parlist.get("Step Acceptance Parameter", 0.05);
109  eta1_ = parlist.get("Radius Shrinking Threshold", 0.05);
110  eta2_ = parlist.get("Radius Growing Threshold", 0.9);
111  gamma0_ = parlist.get("Radius Shrinking Rate (Negative rho)", 0.0625);
112  gamma1_ = parlist.get("Radius Shrinking Rate (Positive rho)", 0.25);
113  gamma2_ = parlist.get("Radius Growing Rate", 2.5);
114  TRsafe_ = parlist.get("Trust-Region Safeguard", 100.0);
116 
117  // Inexactness Information
118  useInexact_.clear();
119  useInexact_.push_back(parlist.get("Use Inexact Objective Function", false));
120  useInexact_.push_back(parlist.get("Use Inexact Gradient", false));
121  useInexact_.push_back(parlist.get("Use Inexact Hessian-Times-A-Vector", false));
122  scale_ = parlist.get("Value Update Tolerance Scaling",1.e-1);
123  omega_ = parlist.get("Value Update Exponent",0.9);
124  force_ = parlist.get("Value Update Forcing Sequence Initial Value",1.0);
125  updateIter_ = parlist.get("Value Update Forcing Sequence Update Frequency",10);
126  forceFactor_ = parlist.get("Value Update Forcing Sequence Reduction Factor",0.1);
127 
128  // Changing Objective Functions
129  softUp_ = parlist.get("Variable Objective Function",false);
130  }
131 
132  virtual void initialize( const Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g) {
133  xupdate_ = x.clone();
134  Hs_ = g.clone();
135  }
136 
137  virtual void update( Vector<Real> &x, Real &fnew, Real &del,
138  int &nfval, int &ngrad, int &flagTR,
139  const Vector<Real> &s, const Real snorm,
140  const Real fold, const Vector<Real> &g,
141  int iter, ProjectedObjective<Real> &pObj ) {
142  Real tol = std::sqrt(ROL_EPSILON);
143 
144  // Compute New Function Value
145  xupdate_->set(x);
146  xupdate_->axpy(1.0,s);
147  /***************************************************************************************************/
148  // BEGIN INEXACT OBJECTIVE FUNCTION COMPUTATION
149  /***************************************************************************************************/
150  Real fold1 = fold;
151  if ( useInexact_[0] ) {
152  if ( !(cnt_%updateIter_) && (cnt_ != 0) ) {
153  force_ *= forceFactor_;
154  }
155  Real c = scale_*std::max(1.e-2,std::min(1.0,1.e4*std::max(pRed_,std::sqrt(ROL_EPSILON))));
156  Real ftol = c*std::pow(std::min(eta1_,1.0-eta2_)
157  *std::min(std::max(pRed_,std::sqrt(ROL_EPSILON)),force_),1.0/omega_);
158  if ( ftol_old_ > ftol || cnt_ == 0 ) {
159  ftol_old_ = ftol;
160  fold1 = pObj.value(x,ftol_old_);
161  }
162  updateObj(*xupdate_,iter,pObj);
163  //pObj.update(*xupdate_,true,iter);
164  fnew = pObj.value(*xupdate_,ftol);
165  cnt_++;
166  }
167  else {
168  updateObj(*xupdate_,iter,pObj);
169  //pObj.update(*xupdate_,true,iter);
170  fnew = pObj.value(*xupdate_,tol);
171  }
172  nfval = 1;
173  Real aRed = fold1 - fnew;
174  /***************************************************************************************************/
175  // FINISH INEXACT OBJECTIVE FUNCTION COMPUTATION
176  /***************************************************************************************************/
177  // If constraints are turned on, then compute a different predicted reduction
178  if (pObj.isConActivated()) {
179  xupdate_->set(x);
180  xupdate_->axpy(-1.0,g.dual());
181  pObj.project(*xupdate_);
182  xupdate_->axpy(-1.0,x);
183  xupdate_->scale(-1.0);
184 
185  pObj.reducedHessVec(*Hs_,s,x,xupdate_->dual(),x,tol);
186  pRed_ = -0.5*s.dot(Hs_->dual());
187 
188  Hs_->set(g);
189  pObj.pruneActive(*Hs_,xupdate_->dual(),x);
190  pRed_ -= s.dot(Hs_->dual());
191  }
192 
193  // Compute Ratio of Actual and Predicted Reduction
194  aRed -= eps_*((1.0 < std::abs(fold1)) ? 1.0 : std::abs(fold1));
195  pRed_ -= eps_*((1.0 < std::abs(fold1)) ? 1.0 : std::abs(fold1));
196  Real rho = 0.0;
197  if ((std::abs(aRed) < eps_) && (std::abs(pRed_) < eps_)) {
198  rho = 1.0;
199  flagTR = 0;
200  }
201  else {
202  rho = aRed/pRed_;
203  if (pRed_ < 0 && aRed > 0) {
204  flagTR = 1;
205  }
206  else if (aRed <= 0 && pRed_ > 0) {
207  flagTR = 2;
208  }
209  else if (aRed <= 0 && pRed_ < 0) {
210  flagTR = 3;
211  }
212  else {
213  flagTR = 0;
214  }
215  }
216 
217  // Check Sufficient Decrease in the Reduced Quadratic Model
218  bool decr = true;
219  if ( pObj.isConActivated() && (std::abs(aRed) > eps_) ) {
220  // Compute Criticality Measure || x - P( x - g ) ||
221  xupdate_->set(x);
222  xupdate_->axpy(-1.0,g.dual());
223  pObj.project(*xupdate_);
224  xupdate_->scale(-1.0);
225  xupdate_->plus(x);
226  Real pgnorm = xupdate_->norm();
227  // Compute Scaled Measure || x - P( x - lam * PI(g) ) ||
228  xupdate_->set(g.dual());
229  pObj.pruneActive(*xupdate_,g,x);
230  Real lam = std::min(1.0, del/xupdate_->norm());
231  xupdate_->scale(-lam);
232  xupdate_->plus(x);
233  pObj.project(*xupdate_);
234  xupdate_->scale(-1.0);
235  xupdate_->plus(x);
236  pgnorm *= xupdate_->norm();
237  // Sufficient decrease?
238  decr = ( aRed >= 0.1*eta0_*pgnorm );
239  flagTR = (!decr ? 4 : flagTR);
240  }
241 
242  // Accept or Reject Step and Update Trust Region
243  if ((rho < eta0_ && flagTR == 0) || flagTR >= 2 || !decr ) { // Step Rejected
244  updateObj(x,iter,pObj);
245  //pObj.update(x,true,iter);
246  fnew = fold1;
247  if (rho < 0.0) { // Negative reduction, interpolate to find new trust-region radius
248  Real gs = s.dot(g.dual());
249  pObj.hessVec(*Hs_,s,x,tol);
250  Real modelVal = s.dot(Hs_->dual());
251  modelVal *= 0.5;
252  modelVal += gs + fold1;
253  Real theta = (1.0-eta2_)*gs/((1.0-eta2_)*(fold1+gs)+eta2_*modelVal-fnew);
254  del = std::min(gamma1_*snorm,std::max(gamma0_,theta)*del);
255  }
256  else { // Shrink trust-region radius
257  del = gamma1_*snorm;
258  }
259  }
260  else if ((rho >= eta0_ && flagTR != 3) || flagTR == 1) { // Step Accepted
261  x.axpy(1.0,s);
262  if (rho >= eta2_) { // Increase trust-region radius
263  del = std::min(gamma2_*del,delmax_);
264  }
265  }
266  }
267 
268  void setPredictedReduction(const Real pRed) {
269  pRed_ = pRed;
270  }
271  Real getPredictedReduction(void) const {
272  return pRed_;
273  }
274 
275  virtual void run( Vector<Real> &s, Real &snorm, Real &del, int &iflag, int &iter, const Vector<Real> &x,
276  const Vector<Real> &grad, const Real &gnorm, ProjectedObjective<Real> &pObj ) = 0;
277 #if 0
278  // If constraints are active, then determine a feasible step
279  if ( pObj.isConActivated() && etr_ != TRUSTREGION_CAUCHYPOINT ) {
280  // Compute projected step stmp = P( x + alpha*s ) - x and xtmp = x + stmp
281  Real alpha = 1.0;
282  Teuchos::RCP<Vector<Real> > stmp = x.clone();
283  stmp->set(s);
284  stmp->scale(alpha);
285  pObj.computeProjectedStep(*stmp,x);
286  Teuchos::RCP<Vector<Real> > xtmp = x.clone();
287  xtmp->set(x);
288  xtmp->axpy(1.0,*stmp);
289  // Compute model components for alpha = 1.0
290  Real tol = std::sqrt(ROL_EPSILON);
291  Teuchos::RCP<Vector<Real> > Bs = x.clone();
292  pObj.hessVec(*Bs,*stmp,x,tol);
293  Real sBs = Bs->dot(*stmp);
294  Real gs = grad.dot(*stmp);
295  Real val = gs + 0.5*sBs;
296  Real val0 = val;
297  // Backtrack alpha until x+alpha*s is feasible
298  int cnt = 0;
299  int maxit = 10;
300  while ( val > val0 || !pObj.isFeasible(*xtmp) ) {
301  // Backtrack alpha
302  alpha *= 0.5;
303  // Computed projected step P( x + alpha*s ) - x and xtmp = x + stmp
304  stmp->set(s);
305  stmp->scale(alpha);
306  pObj.computeProjectedStep(*stmp,x);
307  xtmp->set(x);
308  xtmp->axpy(1.0,*stmp);
309  // Evaluate Model
310  val0 = val;
311  pObj.hessVec(*Bs,*stmp,x,tol);
312  sBs = Bs->dot(*stmp);
313  gs = grad.dot(*stmp);
314  val = gs + 0.5*sBs;
315  // Update iteration count
316  cnt++;
317  if ( cnt >= maxit ) { break; }
318  }
319  s.set(*stmp);
320  pRed_ = -val;
321  }
322 #endif
323 };
324 
325 }
326 
327 #include "ROL_CauchyPoint.hpp"
328 #include "ROL_DogLeg.hpp"
329 #include "ROL_DoubleDogLeg.hpp"
330 #include "ROL_TruncatedCG.hpp"
331 
332 #endif
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:211
Real value(const Vector< Real > &x, Real &tol)
Compute value.
TrustRegion(Teuchos::ParameterList &parlist)
Teuchos::RCP< Vector< Real > > xupdate_
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:141
virtual void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
void reducedHessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &p, const Vector< Real > &d, const Vector< Real > &x, Real &tol)
Apply the reduced Hessian to a vector, v. The reduced Hessian first removes elements of v correspondi...
Contains definitions of custom data types in ROL.
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Provides interface for and implements trust-region subproblem solvers.
void pruneActive(Vector< Real > &v, const Vector< Real > &g, const Vector< Real > &x)
Contains definitions for helper functions in ROL.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
void hessVec(Vector< Real > &Hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual Real dot(const Vector &x) const =0
Compute where .
virtual void run(Vector< Real > &s, Real &snorm, Real &del, int &iflag, int &iter, const Vector< Real > &x, const Vector< Real > &grad, const Real &gnorm, ProjectedObjective< Real > &pObj)=0
Real getPredictedReduction(void) const
Teuchos::RCP< Vector< Real > > Hs_
void updateObj(Vector< Real > &x, int iter, ProjectedObjective< Real > &pObj)
virtual void update(Vector< Real > &x, Real &fnew, Real &del, int &nfval, int &ngrad, int &flagTR, const Vector< Real > &s, const Real snorm, const Real fold, const Vector< Real > &g, int iter, ProjectedObjective< Real > &pObj)
void setPredictedReduction(const Real pRed)
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
void project(Vector< Real > &x)
static const double ROL_OVERFLOW
Platform-dependent maximum double.
Definition: ROL_Types.hpp:123
std::vector< bool > useInexact_
static const double ROL_EPSILON
Platform-dependent machine epsilon.
Definition: ROL_Types.hpp:115