ROL
ROL_TypeU_TrustRegionAlgorithm_Def.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_TRUSTREGIONALGORITHM_U_DEF_H
45 #define ROL_TRUSTREGIONALGORITHM_U_DEF_H
46 
48 
49 namespace ROL {
50 namespace TypeU {
51 
52 template<typename Real>
54  const Ptr<Secant<Real>> &secant )
55  : Algorithm<Real>(), esec_(SECANT_USERDEFINED) {
56  // Set status test
57  status_->reset();
58  status_->add(makePtr<StatusTest<Real>>(parlist));
59 
60  // Trust-Region Parameters
61  ParameterList &slist = parlist.sublist("Step");
62  ParameterList &trlist = slist.sublist("Trust Region");
63  state_->searchSize = trlist.get("Initial Radius", static_cast<Real>(-1));
64  delMax_ = trlist.get("Maximum Radius", ROL_INF<Real>());
65  eta0_ = trlist.get("Step Acceptance Threshold", static_cast<Real>(0.05));
66  eta1_ = trlist.get("Radius Shrinking Threshold", static_cast<Real>(0.05));
67  eta2_ = trlist.get("Radius Growing Threshold", static_cast<Real>(0.9));
68  gamma0_ = trlist.get("Radius Shrinking Rate (Negative rho)", static_cast<Real>(0.0625));
69  gamma1_ = trlist.get("Radius Shrinking Rate (Positive rho)", static_cast<Real>(0.25));
70  gamma2_ = trlist.get("Radius Growing Rate", static_cast<Real>(2.5));
71  TRsafe_ = trlist.get("Safeguard Size", static_cast<Real>(100.0));
72  eps_ = TRsafe_*ROL_EPSILON<Real>();
73  // Nonmonotone Information
74  NMstorage_ = trlist.get("Nonmonotone Storage Limit", 0);
75  useNM_ = (NMstorage_ <= 0 ? false : true);
76  // Inexactness Information
77  ParameterList &glist = parlist.sublist("General");
78  useInexact_.clear();
79  useInexact_.push_back(glist.get("Inexact Objective Function", false));
80  useInexact_.push_back(glist.get("Inexact Gradient", false));
81  useInexact_.push_back(glist.get("Inexact Hessian-Times-A-Vector", false));
82  // Trust-Region Inexactness Parameters
83  ParameterList &ilist = trlist.sublist("Inexact").sublist("Gradient");
84  scale0_ = ilist.get("Tolerance Scaling", static_cast<Real>(0.1));
85  scale1_ = ilist.get("Relative Tolerance", static_cast<Real>(2));
86  // Inexact Function Evaluation Information
87  ParameterList &vlist = trlist.sublist("Inexact").sublist("Value");
88  scale_ = vlist.get("Tolerance Scaling", static_cast<Real>(1.e-1));
89  omega_ = vlist.get("Exponent", static_cast<Real>(0.9));
90  force_ = vlist.get("Forcing Sequence Initial Value", static_cast<Real>(1.0));
91  updateIter_ = vlist.get("Forcing Sequence Update Frequency", static_cast<int>(10));
92  forceFactor_ = vlist.get("Forcing Sequence Reduction Factor", static_cast<Real>(0.1));
93  // Initialize Trust Region Subproblem Solver Object
94  etr_ = StringToETrustRegionU(trlist.get("Subproblem Solver", "Dogleg"));
95  solver_ = TrustRegionUFactory<Real>(parlist);
96  verbosity_ = glist.get("Output Level", 0);
97  // Secant Information
98  useSecantPrecond_ = glist.sublist("Secant").get("Use as Preconditioner", false);
99  useSecantHessVec_ = glist.sublist("Secant").get("Use as Hessian", false);
100  if (secant == nullPtr) {
101  esec_ = StringToESecant(glist.sublist("Secant").get("Type","Limited-Memory BFGS"));
102  }
103  // Initialize trust region model
104  model_ = makePtr<TrustRegionModel_U<Real>>(parlist,secant);
105  printHeader_ = verbosity_ > 2;
106 }
107 
108 template<typename Real>
110  const Vector<Real> &g,
111  Vector<Real> &Bg,
112  Objective<Real> &obj,
113  std::ostream &outStream) {
114  // Initialize data
116  solver_->initialize(x,g);
117  model_->initialize(x,g);
118  // Update approximate gradient and approximate objective function.
119  Real ftol = static_cast<Real>(0.1)*ROL_OVERFLOW<Real>();
120  obj.update(x,UpdateType::Initial,state_->iter);
121  state_->value = obj.value(x,ftol);
122  state_->nfval++;
123  state_->snorm = ROL_INF<Real>();
124  state_->gnorm = ROL_INF<Real>();
125  computeGradient(x,obj);
126  // Check if inverse Hessian is implemented for dogleg methods
127  model_->validate(obj,x,g,etr_);
128  // Compute initial trust region radius if desired.
129  if ( state_->searchSize <= static_cast<Real>(0) ) {
130  int nfval = 0;
131  state_->searchSize
132  = TRUtils::initialRadius<Real>(nfval,x,*state_->gradientVec,Bg,
133  state_->value,state_->gnorm,obj,*model_,delMax_,
134  outStream,(verbosity_>1));
135  state_->nfval += nfval;
136  }
137 }
138 
139 template<typename Real>
141  Objective<Real> &obj,
142  Real pRed) {
143  const Real one(1);
144  Real tol(std::sqrt(ROL_EPSILON<Real>())), fval(0);
145  if ( useInexact_[0] ) {
146  if ( !(state_->iter%updateIter_) && (state_->iter != 0) ) {
147  force_ *= forceFactor_;
148  }
149  Real eta = static_cast<Real>(0.999)*std::min(eta1_,one-eta2_);
150  tol = scale_*std::pow(eta*std::min(pRed,force_),one/omega_);
151  state_->value = obj.value(*state_->iterateVec,tol);
152  state_->nfval++;
153  }
154  // Evaluate objective function at new iterate
155  obj.update(x,UpdateType::Trial);
156  fval = obj.value(x,tol);
157  state_->nfval++;
158  return fval;
159 }
160 
161 template<typename Real>
163  Objective<Real> &obj) {
164  if ( useInexact_[1] ) {
165  const Real one(1);
166  Real gtol1 = scale0_*state_->searchSize;
167  Real gtol0 = gtol1 + one;
168  while ( gtol0 > gtol1 ) {
169  obj.gradient(*state_->gradientVec,x,gtol1);
170  state_->gnorm = state_->gradientVec->norm();
171  gtol0 = gtol1;
172  gtol1 = scale0_*std::min(state_->gnorm,state_->searchSize);
173  }
174  }
175  else {
176  Real gtol = std::sqrt(ROL_EPSILON<Real>());
177  obj.gradient(*state_->gradientVec,x,gtol);
178  state_->gnorm = state_->gradientVec->norm();
179  }
180  state_->ngrad++;
181 }
182 
183 template<typename Real>
185  const Vector<Real> &g,
186  Objective<Real> &obj,
187  std::ostream &outStream ) {
188  const Real zero(0);
189  // Initialize trust-region data
190  Real ftrial(0), pRed(0), rho(0);
191  Ptr<Vector<Real>> gvec = g.clone();
192  initialize(x,g,*gvec,obj,outStream);
193  // Initialize nonmonotone data
194  Real rhoNM(0), sigmac(0), sigmar(0);
195  Real fr(state_->value), fc(state_->value), fmin(state_->value);
196  TRUtils::ETRFlag TRflagNM;
197  int L(0);
198 
199  // Output
200  if (verbosity_ > 0) writeOutput(outStream,true);
201 
202  while (status_->check(*state_)) {
203  // Build trust-region model
204  model_->setData(obj,x,*state_->gradientVec);
205  // Minimize trust-region model over trust-region constraint
206  pRed = zero;
207  SPflag_ = 0; SPiter_ = 0;
208  solver_->solve(*state_->stepVec,state_->snorm,pRed,SPflag_,SPiter_,
209  state_->searchSize,*model_);
210  // Compute trial objective function value
211  x.plus(*state_->stepVec);
212  ftrial = computeValue(x,obj,pRed);
213  // Compute ratio of actual and predicted reduction
214  TRflag_ = TRUtils::SUCCESS;
215  TRUtils::analyzeRatio<Real>(rho,TRflag_,state_->value,ftrial,pRed,eps_,outStream,verbosity_>1);
216  if (useNM_) {
217  TRUtils::analyzeRatio<Real>(rhoNM,TRflagNM,fr,ftrial,pRed+sigmar,eps_,outStream,verbosity_>1);
218  TRflag_ = (rho < rhoNM ? TRflagNM : TRflag_);
219  rho = (rho < rhoNM ? rhoNM : rho );
220  }
221  // Update algorithm state
222  state_->iter++;
223  // Accept/reject step and update trust region radius
224  if ((rho < eta0_ && TRflag_ == TRUtils::SUCCESS)
225  || (TRflag_ >= 2)) { // Step Rejected
226  x.set(*state_->iterateVec);
227  obj.update(x,UpdateType::Revert,state_->iter);
228  if (rho < zero && TRflag_ != TRUtils::TRNAN) {
229  // Negative reduction, interpolate to find new trust-region radius
230  state_->searchSize = TRUtils::interpolateRadius<Real>(*state_->gradientVec,*state_->stepVec,
231  state_->snorm,pRed,state_->value,ftrial,state_->searchSize,gamma0_,gamma1_,eta2_,
232  outStream,verbosity_>1);
233  }
234  else { // Shrink trust-region radius
235  state_->searchSize = gamma1_*std::min(state_->snorm,state_->searchSize);
236  }
237  if (useInexact_[1]) computeGradient(x,obj);
238  }
239  else if ((rho >= eta0_ && TRflag_ != TRUtils::NPOSPREDNEG)
240  || (TRflag_ == TRUtils::POSPREDNEG)) { // Step Accepted
241  state_->iterateVec->set(x);
242  state_->value = ftrial;
243  obj.update(x,UpdateType::Accept,state_->iter);
244  if (useNM_) {
245  sigmac += pRed; sigmar += pRed;
246  if (ftrial < fmin) { fmin = ftrial; fc = fmin; sigmac = zero; L = 0; }
247  else {
248  L++;
249  if (ftrial > fc) { fc = ftrial; sigmac = zero; }
250  if (L == NMstorage_) { fr = fc; sigmar = sigmac; }
251  }
252  }
253  // Increase trust-region radius
254  if (rho >= eta2_) state_->searchSize = std::min(gamma2_*state_->searchSize, delMax_);
255  // Compute gradient at new iterate
256  gvec->set(*state_->gradientVec);
257  computeGradient(x,obj);
258  // Update secant information in trust-region model
259  model_->update(x,*state_->stepVec,*gvec,*state_->gradientVec,
260  state_->snorm,state_->iter);
261  }
262  // Update Output
263  if (verbosity_ > 0) writeOutput(outStream,printHeader_);
264  }
265  if (verbosity_ > 0) Algorithm<Real>::writeExitStatus(outStream);
266 }
267 
268 template<typename Real>
269 void TrustRegionAlgorithm<Real>::writeHeader( std::ostream& os ) const {
270  std::stringstream hist;
271  if(verbosity_ > 1) {
272  hist << std::string(114,'-') << std::endl;
273  hist << "Trust-Region status output definitions" << std::endl << std::endl;
274  hist << " iter - Number of iterates (steps taken)" << std::endl;
275  hist << " value - Objective function value" << std::endl;
276  hist << " gnorm - Norm of the gradient" << std::endl;
277  hist << " snorm - Norm of the step (update to optimization vector)" << std::endl;
278  hist << " delta - Trust-Region radius" << std::endl;
279  hist << " #fval - Number of times the objective function was evaluated" << std::endl;
280  hist << " #grad - Number of times the gradient was computed" << std::endl;
281  hist << std::endl;
282  hist << " tr_flag - Trust-Region flag" << std::endl;
283  for( int flag = TRUtils::SUCCESS; flag != TRUtils::UNDEFINED; ++flag ) {
284  hist << " " << NumberToString(flag) << " - "
285  << TRUtils::ETRFlagToString(static_cast<TRUtils::ETRFlag>(flag)) << std::endl;
286  }
287  if( etr_ == TRUSTREGION_U_TRUNCATEDCG ) {
288  hist << std::endl;
289  hist << " iterCG - Number of Truncated CG iterations" << std::endl << std::endl;
290  hist << " flagGC - Trust-Region Truncated CG flag" << std::endl;
291  for( int flag = CG_FLAG_SUCCESS; flag != CG_FLAG_UNDEFINED; ++flag ) {
292  hist << " " << NumberToString(flag) << " - "
293  << ECGFlagToString(static_cast<ECGFlag>(flag)) << std::endl;
294  }
295  }
296  else if( etr_ == TRUSTREGION_U_SPG ) {
297  hist << std::endl;
298  hist << " iterCG - Number of spectral projected gradient iterations" << std::endl << std::endl;
299  hist << " flagGC - Trust-Region spectral projected gradient flag" << std::endl;
300  }
301  hist << std::string(114,'-') << std::endl;
302  }
303  hist << " ";
304  hist << std::setw(6) << std::left << "iter";
305  hist << std::setw(15) << std::left << "value";
306  hist << std::setw(15) << std::left << "gnorm";
307  hist << std::setw(15) << std::left << "snorm";
308  hist << std::setw(15) << std::left << "delta";
309  hist << std::setw(10) << std::left << "#fval";
310  hist << std::setw(10) << std::left << "#grad";
311  hist << std::setw(10) << std::left << "tr_flag";
312  if ( etr_ == TRUSTREGION_U_TRUNCATEDCG ) {
313  hist << std::setw(10) << std::left << "iterCG";
314  hist << std::setw(10) << std::left << "flagCG";
315  }
316  else if (etr_ == TRUSTREGION_U_SPG) {
317  hist << std::setw(10) << std::left << "iterSPG";
318  hist << std::setw(10) << std::left << "flagSPG";
319  }
320  hist << std::endl;
321  os << hist.str();
322 }
323 
324 template<typename Real>
325 void TrustRegionAlgorithm<Real>::writeName( std::ostream& os ) const {
326  std::stringstream hist;
327  hist << std::endl << ETrustRegionUToString(etr_) << " Trust-Region Solver";
328  if ( useSecantPrecond_ || useSecantHessVec_ ) {
329  if ( useSecantPrecond_ && !useSecantHessVec_ ) {
330  hist << " with " << ESecantToString(esec_) << " Preconditioning" << std::endl;
331  }
332  else if ( !useSecantPrecond_ && useSecantHessVec_ ) {
333  hist << " with " << ESecantToString(esec_) << " Hessian Approximation" << std::endl;
334  }
335  else {
336  hist << " with " << ESecantToString(esec_) << " Preconditioning and Hessian Approximation" << std::endl;
337  }
338  }
339  else {
340  hist << std::endl;
341  }
342  os << hist.str();
343 }
344 
345 template<typename Real>
346 void TrustRegionAlgorithm<Real>::writeOutput(std::ostream& os, bool print_header) const {
347  std::stringstream hist;
348  hist << std::scientific << std::setprecision(6);
349  if ( state_->iter == 0 ) {
350  writeName(os);
351  }
352  if ( print_header ) {
353  writeHeader(os);
354  }
355  if ( state_->iter == 0 ) {
356  hist << " ";
357  hist << std::setw(6) << std::left << state_->iter;
358  hist << std::setw(15) << std::left << state_->value;
359  hist << std::setw(15) << std::left << state_->gnorm;
360  hist << std::setw(15) << std::left << "---";
361  hist << std::setw(15) << std::left << state_->searchSize;
362  hist << std::setw(10) << std::left << state_->nfval;
363  hist << std::setw(10) << std::left << state_->ngrad;
364  hist << std::setw(10) << std::left << "---";
365  if ( etr_ == TRUSTREGION_U_TRUNCATEDCG || etr_ == TRUSTREGION_U_SPG ) {
366  hist << std::setw(10) << std::left << "---";
367  hist << std::setw(10) << std::left << "---";
368  }
369  hist << std::endl;
370  }
371  else {
372  hist << " ";
373  hist << std::setw(6) << std::left << state_->iter;
374  hist << std::setw(15) << std::left << state_->value;
375  hist << std::setw(15) << std::left << state_->gnorm;
376  hist << std::setw(15) << std::left << state_->snorm;
377  hist << std::setw(15) << std::left << state_->searchSize;
378  hist << std::setw(10) << std::left << state_->nfval;
379  hist << std::setw(10) << std::left << state_->ngrad;
380  hist << std::setw(10) << std::left << TRflag_;
381  if ( etr_ == TRUSTREGION_U_TRUNCATEDCG || etr_ == TRUSTREGION_U_SPG ) {
382  hist << std::setw(10) << std::left << SPiter_;
383  hist << std::setw(10) << std::left << SPflag_;
384  }
385  hist << std::endl;
386  }
387  os << hist.str();
388 }
389 } // namespace TypeU
390 } // namespace ROL
391 
392 #endif
std::string ECGFlagToString(ECGFlag cgf)
Definition: ROL_Types.hpp:831
int verbosity_
Print additional information to screen if > 0.
Provides the interface to evaluate objective functions.
void computeGradient(const Vector< Real > &x, Objective< Real > &obj)
Compute gradient to iteratively satisfy inexactness condition.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
void initialize(const Vector< Real > &x, const Vector< Real > &g, Vector< Real > &Bg, Objective< Real > &obj, std::ostream &outStream=std::cout)
virtual void plus(const Vector &x)=0
Compute , where .
const Ptr< AlgorithmState< Real > > state_
Real scale1_
Scale for inexact gradient computation.
ETrustRegionU StringToETrustRegionU(std::string s)
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
TrustRegionAlgorithm(ParameterList &parlist, const Ptr< Secant< Real >> &secant=nullPtr)
void initialize(const Vector< Real > &x, const Vector< Real > &g)
Real delMax_
Maximum trust-region radius.
ESecant StringToESecant(std::string s)
Definition: ROL_Types.hpp:543
Ptr< TrustRegion_U< Real > > solver_
Container for trust-region solver object.
bool printHeader_
Print header at every iteration.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
Contains definitions of enums for trust region algorithms.
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Ptr< TrustRegionModel_U< Real > > model_
Container for trust-region model.
Real gamma0_
Radius decrease rate (negative rho).
ETRFlag
Enumation of flags used by trust-region solvers.
void writeOutput(std::ostream &os, bool print_header=false) const override
Print iterate status.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
std::string NumberToString(T Number)
Definition: ROL_Types.hpp:81
Provides an interface to run unconstrained optimization algorithms.
void run(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, std::ostream &outStream=std::cout) override
Run algorithm on unconstrained problems (Type-U). This general interface supports the use of dual opt...
void writeName(std::ostream &os) const override
Print step name.
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:79
Provides an interface to check status of optimization algorithms.
std::string ETRFlagToString(ETRFlag trf)
virtual void writeExitStatus(std::ostream &os) const
Real scale0_
Scale for inexact gradient computation.
void writeHeader(std::ostream &os) const override
Print iterate header.
Real TRsafe_
Safeguard size for numerically evaluating ratio.
Real gamma1_
Radius decrease rate (positive rho).
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
ETrustRegionU etr_
Trust-region subproblem solver type.
const Ptr< CombinedStatusTest< Real > > status_
std::vector< bool > useInexact_
Flags for inexact (0) objective function, (1) gradient, (2) Hessian.
std::string ESecantToString(ESecant tr)
Definition: ROL_Types.hpp:495
Real eps_
Safeguard for numerically evaluating ratio.
Real computeValue(const Vector< Real > &x, Objective< Real > &obj, Real pRed)
std::string ETrustRegionUToString(ETrustRegionU tr)