ROL
ROL_TrustRegion.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_TRUSTREGION_H
11#define ROL_TRUSTREGION_H
12
16
17#include "ROL_Types.hpp"
22
23namespace ROL {
24
25template<class Real>
27private:
28
29 ROL::Ptr<Vector<Real> > prim_, dual_, xtmp_;
30
32
35 Real pRed_;
37 Real mu0_;
38
39 std::vector<bool> useInexact_;
40
42
45
46 unsigned verbosity_;
47
48 // POST SMOOTHING PARAMETERS
51 Real mu_;
52 Real beta_;
53
54public:
55
56 virtual ~TrustRegion() {}
57
58 // Constructor
59 TrustRegion( ROL::ParameterList &parlist )
60 : pRed_(0), ftol_old_(ROL_OVERFLOW<Real>()), cnt_(0), verbosity_(0) {
61 // Trust-Region Parameters
62 ROL::ParameterList list = parlist.sublist("Step").sublist("Trust Region");
63 std::string modelName = list.get("Subproblem Model", "Kelley-Sachs");
65 eta0_ = list.get("Step Acceptance Threshold", static_cast<Real>(0.05));
66 eta1_ = list.get("Radius Shrinking Threshold", static_cast<Real>(0.05));
67 eta2_ = list.get("Radius Growing Threshold", static_cast<Real>(0.9));
68 gamma0_ = list.get("Radius Shrinking Rate (Negative rho)", static_cast<Real>(0.0625));
69 gamma1_ = list.get("Radius Shrinking Rate (Positive rho)", static_cast<Real>(0.25));
70 gamma2_ = list.get("Radius Growing Rate", static_cast<Real>(2.5));
71 mu0_ = list.get("Sufficient Decrease Parameter", static_cast<Real>(1.e-4));
72 TRsafe_ = list.get("Safeguard Size", static_cast<Real>(100.0));
74 // General Inexactness Information
75 ROL::ParameterList &glist = parlist.sublist("General");
76 useInexact_.clear();
77
78 bool inexactObj = glist.get("Inexact Objective Function", false);
79 bool inexactGrad = glist.get("Inexact Gradient", false);
80 bool inexactHessVec = glist.get("Inexact Hessian-Times-A-Vector", false);
81 useInexact_.push_back(inexactObj );
82 useInexact_.push_back(inexactGrad );
83 useInexact_.push_back(inexactHessVec);
84 // Inexact Function Evaluation Information
85 ROL::ParameterList &ilist = list.sublist("Inexact").sublist("Value");
86 scale_ = ilist.get("Tolerance Scaling", static_cast<Real>(1.e-1));
87 omega_ = ilist.get("Exponent", static_cast<Real>(0.9));
88 force_ = ilist.get("Forcing Sequence Initial Value", static_cast<Real>(1.0));
89 updateIter_ = ilist.get("Forcing Sequence Update Frequency", static_cast<int>(10));
90 forceFactor_ = ilist.get("Forcing Sequence Reduction Factor", static_cast<Real>(0.1));
91 // Get verbosity level
92 verbosity_ = glist.get("Print Verbosity", 0);
93 // Post-smoothing parameters
94 max_fval_ = list.sublist("Post-Smoothing").get("Function Evaluation Limit", 20);
95 alpha_init_ = list.sublist("Post-Smoothing").get("Initial Step Size", static_cast<Real>(1));
96 mu_ = list.sublist("Post-Smoothing").get("Tolerance", static_cast<Real>(0.9999));
97 beta_ = list.sublist("Post-Smoothing").get("Rate", static_cast<Real>(0.01));
98 }
99
100 virtual void initialize( const Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g) {
101 prim_ = x.clone();
102 dual_ = g.clone();
103 xtmp_ = x.clone();
104 }
105
106 virtual void update( Vector<Real> &x,
107 Real &fnew,
108 Real &del,
109 int &nfval,
110 int &ngrad,
111 ETrustRegionFlag &flagTR,
112 const Vector<Real> &s,
113 const Real snorm,
114 const Real fold,
115 const Vector<Real> &g,
116 int iter,
117 Objective<Real> &obj,
119 TrustRegionModel<Real> &model ) {
120 Real tol = std::sqrt(ROL_EPSILON<Real>());
121 const Real one(1), zero(0);
122
123 /***************************************************************************************************/
124 // BEGIN INEXACT OBJECTIVE FUNCTION COMPUTATION
125 /***************************************************************************************************/
126 // Update inexact objective function
127 Real fold1 = fold, ftol = tol; // TOL(1.e-2);
128 if ( useInexact_[0] ) {
129 if ( !(cnt_%updateIter_) && (cnt_ != 0) ) {
131 }
132 //const Real oe4(1e4);
133 //Real c = scale_*std::max(TOL,std::min(one,oe4*std::max(pRed_,std::sqrt(ROL_EPSILON<Real>()))));
134 //ftol = c*std::pow(std::min(eta1_,one-eta2_)
135 // *std::min(std::max(pRed_,std::sqrt(ROL_EPSILON<Real>())),force_),one/omega_);
136 //if ( ftol_old_ > ftol || cnt_ == 0 ) {
137 // ftol_old_ = ftol;
138 // fold1 = obj.value(x,ftol_old_);
139 //}
140 //cnt_++;
141 Real eta = static_cast<Real>(0.999)*std::min(eta1_,one-eta2_);
142 ftol = scale_*std::pow(eta*std::min(pRed_,force_),one/omega_);
143 ftol_old_ = ftol;
144 fold1 = obj.value(x,ftol_old_);
145 cnt_++;
146 }
147 // Evaluate objective function at new iterate
148 prim_->set(x); prim_->plus(s);
149 if (bnd.isActivated()) {
150 bnd.project(*prim_);
151 }
152 obj.update(*prim_);
153 fnew = obj.value(*prim_,ftol);
154
155 nfval = 1;
156 Real aRed = fold1 - fnew;
157 /***************************************************************************************************/
158 // FINISH INEXACT OBJECTIVE FUNCTION COMPUTATION
159 /***************************************************************************************************/
160
161 /***************************************************************************************************/
162 // BEGIN COMPUTE RATIO OF ACTUAL AND PREDICTED REDUCTION
163 /***************************************************************************************************/
164 // Modify Actual and Predicted Reduction According to Model
165 model.updateActualReduction(aRed,s);
167
168 if ( verbosity_ > 0 ) {
169 std::cout << std::endl;
170 std::cout << " Computation of actual and predicted reduction" << std::endl;
171 std::cout << " Current objective function value: " << fold1 << std::endl;
172 std::cout << " New objective function value: " << fnew << std::endl;
173 std::cout << " Actual reduction: " << aRed << std::endl;
174 std::cout << " Predicted reduction: " << pRed_ << std::endl;
175 }
176
177 // Compute Ratio of Actual and Predicted Reduction
178 Real EPS = eps_*((one > std::abs(fold1)) ? one : std::abs(fold1));
179 Real aRed_safe = aRed + EPS, pRed_safe = pRed_ + EPS;
180 Real rho(0);
181 if (((std::abs(aRed_safe) < eps_) && (std::abs(pRed_safe) < eps_)) || aRed == pRed_) {
182 rho = one;
184 }
185 else if ( std::isnan(aRed_safe) || std::isnan(pRed_safe) ) {
186 rho = -one;
187 flagTR = TRUSTREGION_FLAG_NAN;
188 }
189 else {
190 rho = aRed_safe/pRed_safe;
191 if (pRed_safe < zero && aRed_safe > zero) {
193 }
194 else if (aRed_safe <= zero && pRed_safe > zero) {
196 }
197 else if (aRed_safe <= zero && pRed_safe < zero) {
199 }
200 else {
202 }
203 }
204
205 if ( verbosity_ > 0 ) {
206 std::cout << " Safeguard: " << eps_ << std::endl;
207 std::cout << " Actual reduction with safeguard: " << aRed_safe << std::endl;
208 std::cout << " Predicted reduction with safeguard: " << pRed_safe << std::endl;
209 std::cout << " Ratio of actual and predicted reduction: " << rho << std::endl;
210 std::cout << " Trust-region flag: " << flagTR << std::endl;
211 }
212 /***************************************************************************************************/
213 // FINISH COMPUTE RATIO OF ACTUAL AND PREDICTED REDUCTION
214 /***************************************************************************************************/
215
216 /***************************************************************************************************/
217 // BEGIN CHECK SUFFICIENT DECREASE FOR BOUND CONSTRAINED PROBLEMS
218 /***************************************************************************************************/
219 bool decr = true;
221 if ( rho >= eta0_ && (std::abs(aRed_safe) > eps_) ) {
222 // Compute Criticality Measure || x - P( x - g ) ||
223 prim_->set(x);
224 prim_->axpy(-one,g.dual());
225 bnd.project(*prim_);
226 prim_->scale(-one);
227 prim_->plus(x);
228 Real pgnorm = prim_->norm();
229 // Compute Scaled Measure || x - P( x - lam * PI(g) ) ||
230 prim_->set(g.dual());
231 bnd.pruneActive(*prim_,g,x);
232 Real lam = std::min(one, del/prim_->norm());
233 prim_->scale(-lam);
234 prim_->plus(x);
235 bnd.project(*prim_);
236 prim_->scale(-one);
237 prim_->plus(x);
238 pgnorm *= prim_->norm();
239 // Sufficient decrease?
240 decr = ( aRed_safe >= mu0_*pgnorm );
241 flagTR = (!decr ? TRUSTREGION_FLAG_QMINSUFDEC : flagTR);
242
243 if ( verbosity_ > 0 ) {
244 std::cout << " Decrease lower bound (constraints): " << mu0_*pgnorm << std::endl;
245 std::cout << " Trust-region flag (constraints): " << flagTR << std::endl;
246 std::cout << " Is step feasible: " << bnd.isFeasible(x) << std::endl;
247 }
248 }
249 }
250 /***************************************************************************************************/
251 // FINISH CHECK SUFFICIENT DECREASE FOR BOUND CONSTRAINED PROBLEMS
252 /***************************************************************************************************/
253
254 /***************************************************************************************************/
255 // BEGIN STEP ACCEPTANCE AND TRUST REGION RADIUS UPDATE
256 /***************************************************************************************************/
257 if ( verbosity_ > 0 ) {
258 std::cout << " Norm of step: " << snorm << std::endl;
259 std::cout << " Trust-region radius before update: " << del << std::endl;
260 }
261 if ((rho < eta0_ && flagTR == TRUSTREGION_FLAG_SUCCESS) || flagTR >= 2 || !decr ) { // Step Rejected
262 fnew = fold1; // This is a bug if rho < zero...
263 if (rho < zero) { // Negative reduction, interpolate to find new trust-region radius
264 Real gs(0);
265 if ( bnd.isActivated() ) {
266 model.dualTransform(*dual_, *model.getGradient());
267 gs = dual_->dot(s.dual());
268 }
269 else {
270 gs = g.dot(s.dual());
271 }
272 Real modelVal = model.value(s,tol);
273 modelVal += fold1;
274 Real theta = (one-eta2_)*gs/((one-eta2_)*(fold1+gs)+eta2_*modelVal-fnew);
275 del = std::min(gamma1_*std::min(snorm,del),std::max(gamma0_,theta)*del);
276 if ( verbosity_ > 0 ) {
277 std::cout << " Interpolation model value: " << modelVal << std::endl;
278 std::cout << " Interpolation step length: " << theta << std::endl;
279 }
280 }
281 else { // Shrink trust-region radius
282 del = gamma1_*std::min(snorm,del);
283 }
284 obj.update(x,true,iter);
285 }
286 else if ((rho >= eta0_ && flagTR != TRUSTREGION_FLAG_NPOSPREDNEG) ||
287 (flagTR == TRUSTREGION_FLAG_POSPREDNEG)) { // Step Accepted
288 // Perform line search (smoothing) to ensure decrease
290 // Compute new gradient
291 xtmp_->set(x); xtmp_->plus(s);
292 bnd.project(*xtmp_);
293 obj.update(*xtmp_);
294 obj.gradient(*dual_,*xtmp_,tol); // MUST DO SOMETHING HERE WITH TOL
295 ngrad++;
296 // Compute smoothed step
297 Real alpha(1);
298 prim_->set(*xtmp_);
299 prim_->axpy(-alpha/alpha_init_,dual_->dual());
300 bnd.project(*prim_);
301 // Compute new objective value
302 obj.update(*prim_);
303 Real ftmp = obj.value(*prim_,tol); // MUST DO SOMETHING HERE WITH TOL
304 nfval++;
305 // Perform smoothing
306 int cnt = 0;
307 alpha = alpha_init_;
308 while ( (ftmp-fnew) >= mu_*aRed ) {
309 prim_->set(*xtmp_);
310 prim_->axpy(-alpha/alpha_init_,dual_->dual());
311 bnd.project(*prim_);
312 obj.update(*prim_);
313 ftmp = obj.value(*prim_,tol); // MUST DO SOMETHING HERE WITH TOL
314 nfval++;
315 if ( cnt >= max_fval_ ) {
316 break;
317 }
318 alpha *= beta_;
319 cnt++;
320 }
321 // Store objective function and iteration information
322 if (std::isnan(ftmp)) {
323 flagTR = TRUSTREGION_FLAG_NAN;
324 del = gamma1_*std::min(snorm,del);
325 rho = static_cast<Real>(-1);
326 //x.axpy(static_cast<Real>(-1),s);
327 //obj.update(x,true,iter);
328 fnew = fold1;
329 }
330 else {
331 fnew = ftmp;
332 x.set(*prim_);
333 }
334 }
335 else {
336 x.plus(s);
337 }
338 if (rho >= eta2_) { // Increase trust-region radius
339 del = gamma2_*del;
340 }
341 obj.update(x,true,iter);
342 }
343
344 if ( verbosity_ > 0 ) {
345 std::cout << " Trust-region radius after update: " << del << std::endl;
346 std::cout << std::endl;
347 }
348 /***************************************************************************************************/
349 // FINISH STEP ACCEPTANCE AND TRUST REGION RADIUS UPDATE
350 /***************************************************************************************************/
351 }
352
353 virtual void run( Vector<Real> &s, // Step (to be computed)
354 Real &snorm, // Step norm (to be computed)
355 int &iflag, // Exit flag (to be computed)
356 int &iter, // Iteration count (to be computed)
357 const Real del, // Trust-region radius
358 TrustRegionModel<Real> &model ) = 0; // Trust-region model
359
360 void setPredictedReduction(const Real pRed) {
361 pRed_ = pRed;
362 }
363
364 Real getPredictedReduction(void) const {
365 return pRed_;
366 }
367};
368
369}
370
372
373#endif
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0 zero)()
Contains definitions of enums for trust region algorithms.
Contains definitions of custom data types in ROL.
Provides the interface to apply upper and lower bound constraints.
virtual bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the -active set.
bool isActivated(void) const
Check if bounds are on.
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
Provides the interface to evaluate objective functions.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
Provides the interface to evaluate trust-region model functions.
virtual void updatePredictedReduction(Real &pred, const Vector< Real > &s)
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
virtual const Ptr< const Vector< Real > > getGradient(void) const
virtual void updateActualReduction(Real &ared, const Vector< Real > &s)
virtual Real value(const Vector< Real > &s, Real &tol)
ETrustRegionModel TRmodel_
Real mu_
Post-Smoothing tolerance for projected methods.
Real getPredictedReduction(void) const
ROL::Ptr< Vector< Real > > xtmp_
TrustRegion(ROL::ParameterList &parlist)
virtual void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
Real alpha_init_
Initial line-search parameter for projected methods.
ROL::Ptr< Vector< Real > > dual_
ROL::Ptr< Vector< Real > > prim_
virtual void run(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)=0
std::vector< bool > useInexact_
void setPredictedReduction(const Real pRed)
Real beta_
Post-Smoothing rate for projected methods.
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)
int max_fval_
Maximum function evaluations in line-search for projected methods.
Defines the linear algebra or vector space interface.
virtual void set(const Vector &x)
Set where .
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
virtual void plus(const Vector &x)=0
Compute , where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual Real dot(const Vector &x) const =0
Compute where .
ETrustRegionModel StringToETrustRegionModel(std::string s)
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
Definition ROL_Types.hpp:57
Real ROL_OVERFLOW(void)
Platform-dependent maximum double.
Definition ROL_Types.hpp:68
@ TRUSTREGION_FLAG_POSPREDNEG
@ TRUSTREGION_FLAG_NPOSPREDNEG
@ TRUSTREGION_FLAG_QMINSUFDEC
@ TRUSTREGION_FLAG_NPOSPREDPOS
@ TRUSTREGION_MODEL_KELLEYSACHS