ROL
ROL_TrustRegionStep.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_TRUSTREGIONSTEP_H
45 #define ROL_TRUSTREGIONSTEP_H
46 
47 #include "ROL_Step.hpp"
48 #include "ROL_Types.hpp"
49 #include "ROL_Secant.hpp"
50 #include "ROL_TrustRegion.hpp"
51 #include <sstream>
52 #include <iomanip>
53 
126 namespace ROL {
127 
128 template <class Real>
129 class TrustRegionStep : public Step<Real> {
130 private:
131 
132  // ADDITIONAL VECTOR STORAGE
133  Ptr<Vector<Real>> xnew_;
134  Ptr<Vector<Real>> xold_;
135  Ptr<Vector<Real>> gp_;
136 
137  // TRUST REGION INFORMATION
138  Ptr<TrustRegion<Real>> trustRegion_;
139  Ptr<TrustRegionModel<Real>> model_;
142  Real delMax_;
144  int SPflag_;
145  int SPiter_;
146  bool bndActive_;
147 
148  // SECANT INFORMATION
149  Ptr<Secant<Real>> secant_;
153 
154  // BOUND CONSTRAINED PARAMETERS
155  Real scaleEps_;
157 
158  // POST SMOOTHING PARAMETERS
159  Real alpha_init_;
160  int max_fval_;
161  Real mu_;
162  Real beta_;
163 
164  // COLEMAN-LI PARAMETERS
168 
169  // INEXACT COMPUTATION PARAMETERS
170  std::vector<bool> useInexact_;
171  Real scale0_;
172  Real scale1_;
173 
174  // VERBOSITY SETTING
176 
183  void parseParameterList(ROL::ParameterList &parlist) {
184  ROL::Ptr<StepState<Real>> step_state = Step<Real>::getState();
185  // Trust-Region Parameters
186  ROL::ParameterList &slist = parlist.sublist("Step");
187  ROL::ParameterList &list = slist.sublist("Trust Region");
188  step_state->searchSize = list.get("Initial Radius", static_cast<Real>(-1));
189  delMax_ = list.get("Maximum Radius", static_cast<Real>(1.e8));
190  // Inexactness Information
191  ROL::ParameterList &glist = parlist.sublist("General");
192  useInexact_.clear();
193  useInexact_.push_back(glist.get("Inexact Objective Function", false));
194  useInexact_.push_back(glist.get("Inexact Gradient", false));
195  useInexact_.push_back(glist.get("Inexact Hessian-Times-A-Vector", false));
196  // Trust-Region Inexactness Parameters
197  ROL::ParameterList &ilist = list.sublist("Inexact").sublist("Gradient");
198  scale0_ = ilist.get("Tolerance Scaling", static_cast<Real>(0.1));
199  scale1_ = ilist.get("Relative Tolerance", static_cast<Real>(2));
200  // Initialize Trust Region Subproblem Solver Object
201  etr_ = StringToETrustRegion(list.get("Subproblem Solver", "Dogleg"));
202  TRmodel_ = StringToETrustRegionModel(list.get("Subproblem Model", "Kelley-Sachs"));
203  useProjectedGrad_ = glist.get("Projected Gradient Criticality Measure", false);
204  trustRegion_ = TrustRegionFactory<Real>(parlist);
205  // Scale for epsilon active sets
206  scaleEps_ = glist.get("Scale for Epsilon Active Sets", static_cast<Real>(1));
207  verbosity_ = glist.get("Print Verbosity", 0);
208  // Post-smoothing parameters
209  max_fval_ = list.sublist("Post-Smoothing").get("Function Evaluation Limit", 20);
210  alpha_init_ = list.sublist("Post-Smoothing").get("Initial Step Size", static_cast<Real>(1));
211  mu_ = list.sublist("Post-Smoothing").get("Tolerance", static_cast<Real>(0.9999));
212  beta_ = list.sublist("Post-Smoothing").get("Rate", static_cast<Real>(0.01));
213  // Coleman-Li parameters
214  stepBackMax_ = list.sublist("Coleman-Li").get("Maximum Step Back", static_cast<Real>(0.9999));
215  stepBackScale_ = list.sublist("Coleman-Li").get("Maximum Step Scale", static_cast<Real>(1));
216  singleReflect_ = list.sublist("Coleman-Li").get("Single Reflection", true);
217  }
218 
234  AlgorithmState<Real> &algo_state ) {
235  Ptr<StepState<Real>> state = Step<Real>::getState();
236  if ( useInexact_[1] ) {
237  const Real one(1);
238  //const Real oem2(1.e-2), oe4(1.e4);
239  //Real c = scale0_*std::max(oem2,std::min(one,oe4*algo_state.gnorm));
240  //Real gtol1 = c*std::min(algo_state.gnorm,state->searchSize);
241  //Real gtol0 = scale1_*gtol1 + one;
242  //while ( gtol0 > gtol1*scale1_ ) {
243  // obj.gradient(*(state->gradientVec),x,gtol1);
244  // algo_state.gnorm = computeCriticalityMeasure(*(state->gradientVec),x,bnd);
245  // gtol0 = gtol1;
246  // c = scale0_*std::max(oem2,std::min(one,oe4*algo_state.gnorm));
247  // gtol1 = c*std::min(algo_state.gnorm,state->searchSize);
248  //}
249  //algo_state.ngrad++;
250  Real gtol1 = scale0_*state->searchSize;
251  //Real gtol1 = scale0_*std::min(algo_state.gnorm,state->searchSize);
252  Real gtol0 = gtol1 + one;
253  while ( gtol0 > gtol1 ) {
254  obj.gradient(*(state->gradientVec),x,gtol1);
255  algo_state.gnorm = computeCriticalityMeasure(*(state->gradientVec),x,bnd);
256  gtol0 = gtol1;
257  gtol1 = scale0_*std::min(algo_state.gnorm,state->searchSize);
258  }
259  algo_state.ngrad++;
260  }
261  else {
262  Real gtol = std::sqrt(ROL_EPSILON<Real>());
263  obj.gradient(*(state->gradientVec),x,gtol);
264  algo_state.ngrad++;
265  algo_state.gnorm = computeCriticalityMeasure(*(state->gradientVec),x,bnd);
266  }
267  }
268 
278  if ( bnd.isActivated() ) {
279  if ( useProjectedGrad_ ) {
280  gp_->set(g);
281  bnd.computeProjectedGradient( *gp_, x );
282  return gp_->norm();
283  }
284  else {
285  Real one(1);
286  xnew_->set(x);
287  xnew_->axpy(-one,g.dual());
288  bnd.project(*xnew_);
289  xnew_->axpy(-one,x);
290  return xnew_->norm();
291  }
292  }
293  else {
294  return g.norm();
295  }
296  }
297 
298 public:
299 
301  using Step<Real>::compute;
302  using Step<Real>::update;
303 
304  virtual ~TrustRegionStep() {}
305 
313  TrustRegionStep( ROL::ParameterList & parlist )
314  : Step<Real>(),
315  xnew_(nullPtr), xold_(nullPtr), gp_(nullPtr),
316  trustRegion_(nullPtr), model_(nullPtr),
319  SPflag_(0), SPiter_(0), bndActive_(false),
320  secant_(nullPtr), esec_(SECANT_LBFGS),
321  useSecantHessVec_(false), useSecantPrecond_(false),
322  scaleEps_(1), useProjectedGrad_(false),
323  alpha_init_(1), max_fval_(20), mu_(0.9999), beta_(0.01),
324  stepBackMax_(0.9999), stepBackScale_(1), singleReflect_(true),
325  scale0_(1), scale1_(1),
326  verbosity_(0) {
327  // Parse input parameterlist
328  parseParameterList(parlist);
329  // Create secant object
330  ROL::ParameterList &glist = parlist.sublist("General");
331  esec_ = StringToESecant(glist.sublist("Secant").get("Type","Limited-Memory BFGS"));
332  useSecantPrecond_ = glist.sublist("Secant").get("Use as Preconditioner", false);
333  useSecantHessVec_ = glist.sublist("Secant").get("Use as Hessian", false);
334  secant_ = SecantFactory<Real>(parlist);
335  }
336 
346  TrustRegionStep( ROL::Ptr<Secant<Real> > &secant, ROL::ParameterList &parlist )
347  : Step<Real>(),
348  xnew_(nullPtr), xold_(nullPtr), gp_(nullPtr),
349  trustRegion_(nullPtr), model_(nullPtr),
352  SPflag_(0), SPiter_(0), bndActive_(false),
353  secant_(nullPtr), esec_(SECANT_LBFGS),
354  useSecantHessVec_(false), useSecantPrecond_(false),
355  scaleEps_(1), useProjectedGrad_(false),
356  alpha_init_(1), max_fval_(20), mu_(0.9999), beta_(0.01),
357  stepBackMax_(0.9999), stepBackScale_(1), singleReflect_(true),
358  scale0_(1), scale1_(1),
359  verbosity_(0) {
360  // Parse input parameterlist
361  parseParameterList(parlist);
362  // Create secant object
363  ROL::ParameterList &glist = parlist.sublist("General");
364  useSecantPrecond_ = glist.sublist("Secant").get("Use as Preconditioner", false);
365  useSecantHessVec_ = glist.sublist("Secant").get("Use as Hessian", false);
366  if ( ROL::is_nullPtr(secant_) ) {
367  ROL::ParameterList Slist;
368  Slist.sublist("General").sublist("Secant").set("Type","Limited-Memory BFGS");
369  Slist.sublist("General").sublist("Secant").set("Maximum Storage",10);
370  secant_ = SecantFactory<Real>(Slist);
371  }
372  }
373 
382  void initialize( Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g,
384  AlgorithmState<Real> &algo_state ) {
386  throw Exception::NotImplemented(">>> ROL::TrustRegionStep : Invalid Trust Region Solver and Model pair!");
387  }
388  Real p1(0.1), oe10(1.e10), zero(0), one(1), half(0.5), three(3), two(2), six(6);
389  Ptr<StepState<Real>> step_state = Step<Real>::getState();
390  bndActive_ = bnd.isActivated();
391 
392  trustRegion_->initialize(x,s,g);
393 
394  Real htol = std::sqrt(ROL_EPSILON<Real>());
395  Real ftol = p1*ROL_OVERFLOW<Real>();
396 
397  step_state->descentVec = s.clone();
398  step_state->gradientVec = g.clone();
399 
400  if ( bnd.isActivated() ) {
401  // Make initial guess feasible
403  bnd.projectInterior(x);
404  }
405  else {
406  bnd.project(x);
407  }
408  xnew_ = x.clone();
409  xold_ = x.clone();
410  }
411  gp_ = g.clone();
412 
413  // Update approximate gradient and approximate objective function.
414  obj.update(x,true,algo_state.iter);
415  algo_state.snorm = oe10;
416  algo_state.value = obj.value(x,ftol);
417  algo_state.nfval++;
418  algo_state.gnorm = ROL_INF<Real>();
419  updateGradient(x,obj,bnd,algo_state);
420 
421  // Try to apply inverse Hessian
422  if ( !useSecantHessVec_ &&
424  try {
425  Ptr<Vector<Real>> v = g.clone();
426  Ptr<Vector<Real>> hv = x.clone();
427  obj.invHessVec(*hv,*v,x,htol);
428  }
429  catch (std::exception &e) {
430  useSecantHessVec_ = true;
431  }
432  }
433 
434  // Evaluate Objective Function at Cauchy Point
435  bool autoRad = false;
436  if ( step_state->searchSize <= zero ) {
437  autoRad = true;
438  Ptr<Vector<Real>> Bg = g.clone();
439  if ( useSecantHessVec_ ) {
440  secant_->applyB(*Bg,(step_state->gradientVec)->dual());
441  }
442  else {
443  obj.hessVec(*Bg,(step_state->gradientVec)->dual(),x,htol);
444  }
445  Real gBg = Bg->dot(*(step_state->gradientVec));
446  Real alpha = one;
447  if ( gBg > ROL_EPSILON<Real>() ) {
448  alpha = algo_state.gnorm*algo_state.gnorm/gBg;
449  }
450  // Evaluate the objective function at the Cauchy point
451  Ptr<Vector<Real>> cp = s.clone();
452  cp->set((step_state->gradientVec)->dual());
453  cp->scale(-alpha);
454  Ptr<Vector<Real>> xcp = x.clone();
455  xcp->set(x);
456  xcp->plus(*cp);
457  if ( bnd.isActivated() ) {
458  bnd.project(*xcp);
459  }
460  obj.update(*xcp);
461  Real fnew = obj.value(*xcp,ftol); // MUST DO SOMETHING HERE WITH FTOL
462  algo_state.nfval++;
463  // Perform cubic interpolation to determine initial trust region radius
464  Real gs = cp->dot((step_state->gradientVec)->dual());
465  Real a = fnew - algo_state.value - gs - half*alpha*alpha*gBg;
466  if ( std::abs(a) < ROL_EPSILON<Real>() ) {
467  // a = 0 implies the objective is quadratic in the negative gradient direction
468  step_state->searchSize = std::min(alpha*algo_state.gnorm,delMax_);
469  }
470  else {
471  Real b = half*alpha*alpha*gBg;
472  Real c = gs;
473  if ( b*b-three*a*c > ROL_EPSILON<Real>() ) {
474  // There is at least one critical point
475  Real t1 = (-b-std::sqrt(b*b-three*a*c))/(three*a);
476  Real t2 = (-b+std::sqrt(b*b-three*a*c))/(three*a);
477  if ( six*a*t1 + two*b > zero ) {
478  // t1 is the minimizer
479  step_state->searchSize = std::min(t1*alpha*algo_state.gnorm,delMax_);
480  }
481  else {
482  // t2 is the minimizer
483  step_state->searchSize = std::min(t2*alpha*algo_state.gnorm,delMax_);
484  }
485  }
486  else {
487  step_state->searchSize = std::min(alpha*algo_state.gnorm,delMax_);
488  }
489  }
490  if (step_state->searchSize <= ROL_EPSILON<Real>()*algo_state.gnorm && autoRad) {
491  step_state->searchSize = one;
492  }
493  }
494  // Build trust-region model
495  if (bnd.isActivated()) {
497  model_ = makePtr<KelleySachsModel<Real>>(obj,
498  bnd,
499  x,
500  *(step_state->gradientVec),
501  secant_,
504  }
505  else if ( TRmodel_ == TRUSTREGION_MODEL_COLEMANLI ) {
506  model_ = makePtr<ColemanLiModel<Real>>(obj,
507  bnd,
508  x,
509  *(step_state->gradientVec),
510  stepBackMax_,
513  secant_,
516  }
517  else if ( TRmodel_ == TRUSTREGION_MODEL_LINMORE ) {
518  model_ = makePtr<LinMoreModel<Real>>(obj,
519  bnd,
520  x,
521  *(step_state->gradientVec),
522  secant_,
525  }
526  else {
527  ROL_TEST_FOR_EXCEPTION( true, std::invalid_argument,
528  ">>> ERROR (TrustRegionStep): Invalid trust-region model!");
529  }
530  }
531  else {
532  model_ = makePtr<TrustRegionModel<Real>>(obj,
533  bnd,
534  x,
535  *(step_state->gradientVec),
536  secant_,
539  }
540  }
541 
552  void compute( Vector<Real> &s, const Vector<Real> &x,
554  AlgorithmState<Real> &algo_state ) {
555  // Get step state
556  Ptr<StepState<Real>> step_state = Step<Real>::getState();
557  // Build trust-region model
558  model_->update(obj,bnd,x,*step_state->gradientVec,secant_);
559  if (bnd.isActivated()) {
561 // Real eps = scaleEps_*algo_state.gnorm;
562  Real eps = scaleEps_ * std::min(std::pow(algo_state.gnorm,static_cast<Real>(0.75)),
563  static_cast<Real>(0.001));
564  dynamicPtrCast<KelleySachsModel<Real>>(model_)->setEpsilon(eps);
565  }
566  else if ( TRmodel_ == TRUSTREGION_MODEL_COLEMANLI ) {
567  dynamicPtrCast<ColemanLiModel<Real>>(model_)->setRadius(step_state->searchSize);
568  }
569  }
570  // Minimize trust-region model over trust-region constraint
571  SPflag_ = 0; SPiter_ = 0;
572  trustRegion_->run(s,algo_state.snorm,SPflag_,SPiter_,step_state->searchSize,*model_);
573  }
574 
587  const Vector<Real> &s,
588  Objective<Real> &obj,
590  AlgorithmState<Real> &algo_state ) {
591  // Get step state
592  Ptr<StepState<Real>> state = Step<Real>::getState();
593  // Store previous step for constraint computations
594  if ( bnd.isActivated() ) {
595  xold_->set(x);
596  }
597  // Update trust-region information;
598  // Performs a hard update on the objective function
600  state->nfval = 0;
601  state->ngrad = 0;
602  Real fold = algo_state.value;
603  Real fnew(0);
604  algo_state.iter++;
605  trustRegion_->update(x,fnew,state->searchSize,state->nfval,state->ngrad,TRflag_,
606  s,algo_state.snorm,fold,*(state->gradientVec),algo_state.iter,
607  obj,bnd,*model_);
608  algo_state.nfval += state->nfval;
609  algo_state.ngrad += state->ngrad;
610  state->flag = static_cast<int>(TRflag_);
611  state->SPiter = SPiter_;
612  state->SPflag = SPflag_;
613  // If step is accepted ...
614  // Compute new gradient and update secant storage
617  // Perform line search (smoothing) to ensure decrease
619  Real tol = std::sqrt(ROL_EPSILON<Real>());
620  // Compute new gradient
621  obj.gradient(*gp_,x,tol); // MUST DO SOMETHING HERE WITH TOL
622  algo_state.ngrad++;
623  // Compute smoothed step
624  Real alpha(1);
625  xnew_->set(x);
626  xnew_->axpy(-alpha*alpha_init_,gp_->dual());
627  bnd.project(*xnew_);
628  // Compute new objective value
629  obj.update(*xnew_,true,algo_state.iter);
630  Real ftmp = obj.value(*xnew_,tol); // MUST DO SOMETHING HERE WITH TOL
631  algo_state.nfval++;
632  // Perform smoothing
633  int cnt = 0;
634  alpha = static_cast<Real>(1)/alpha_init_;
635  while ( (fnew-ftmp) <= mu_*(fnew-fold) ) {
636  xnew_->set(x);
637  xnew_->axpy(-alpha*alpha_init_,gp_->dual());
638  bnd.project(*xnew_);
639  obj.update(*xnew_,true,algo_state.iter);
640  ftmp = obj.value(*xnew_,tol); // MUST DO SOMETHING HERE WITH TOL
641  algo_state.nfval++;
642  if ( cnt >= max_fval_ ) {
643  break;
644  }
645  alpha *= beta_;
646  cnt++;
647  }
648  // Store objective function and iteration information
649  fnew = ftmp;
650  x.set(*xnew_);
651  }
652  // Store previous gradient for secant update
654  gp_->set(*(state->gradientVec));
655  }
656  // Update objective function and approximate model
657  updateGradient(x,obj,bnd,algo_state);
658  // Update secant information
660  if ( bnd.isActivated() ) { // Compute new constrained step
661  xnew_->set(x);
662  xnew_->axpy(-static_cast<Real>(1),*xold_);
663  secant_->updateStorage(x,*(state->gradientVec),*gp_,*xnew_,algo_state.snorm,algo_state.iter+1);
664  }
665  else {
666  secant_->updateStorage(x,*(state->gradientVec),*gp_,s,algo_state.snorm,algo_state.iter+1);
667  }
668  }
669  // Update algorithm state
670  (algo_state.iterateVec)->set(x);
671  }
672  else {
673  if ( useInexact_[1] ) {
674  // Update objective function and approximate model
675  updateGradient(x,obj,bnd,algo_state);
676  }
677  }
678  // Update algorithm state
679  algo_state.value = fnew;
680  }
681 
686  std::string printHeader( void ) const {
687  std::stringstream hist;
688 
689  if(verbosity_>0) {
690  hist << std::string(114,'-') << "\n";
691 
692  hist << "Trust-Region status output definitions\n\n";
693 
694  hist << " iter - Number of iterates (steps taken) \n";
695  hist << " value - Objective function value \n";
696  hist << " gnorm - Norm of the gradient\n";
697  hist << " snorm - Norm of the step (update to optimization vector)\n";
698  hist << " delta - Trust-Region radius\n";
699  hist << " #fval - Number of times the objective function was evaluated\n";
700  hist << " #grad - Number of times the gradient was computed\n";
701 
702 
703 
704  hist << "\n";
705  hist << " tr_flag - Trust-Region flag" << "\n";
706  for( int flag = TRUSTREGION_FLAG_SUCCESS; flag != TRUSTREGION_FLAG_UNDEFINED; ++flag ) {
707  hist << " " << NumberToString(flag) << " - "
708  << ETrustRegionFlagToString(static_cast<ETrustRegionFlag>(flag)) << "\n";
709 
710  }
711 
712  if( etr_ == TRUSTREGION_TRUNCATEDCG ) {
713  hist << "\n";
714  hist << " iterCG - Number of Truncated CG iterations\n\n";
715  hist << " flagGC - Trust-Region Truncated CG flag" << "\n";
716  for( int flag = CG_FLAG_SUCCESS; flag != CG_FLAG_UNDEFINED; ++flag ) {
717  hist << " " << NumberToString(flag) << " - "
718  << ECGFlagToString(static_cast<ECGFlag>(flag)) << "\n";
719  }
720  }
721 
722  hist << std::string(114,'-') << "\n";
723  }
724 
725  hist << " ";
726  hist << std::setw(6) << std::left << "iter";
727  hist << std::setw(15) << std::left << "value";
728  hist << std::setw(15) << std::left << "gnorm";
729  hist << std::setw(15) << std::left << "snorm";
730  hist << std::setw(15) << std::left << "delta";
731  hist << std::setw(10) << std::left << "#fval";
732  hist << std::setw(10) << std::left << "#grad";
733  hist << std::setw(10) << std::left << "tr_flag";
735  hist << std::setw(10) << std::left << "iterCG";
736  hist << std::setw(10) << std::left << "flagCG";
737  }
738  hist << "\n";
739  return hist.str();
740  }
741 
746  std::string printName( void ) const {
747  std::stringstream hist;
748  hist << "\n" << ETrustRegionToString(etr_) << " Trust-Region Solver";
751  hist << " with " << ESecantToString(esec_) << " Preconditioning\n";
752  }
753  else if ( !useSecantPrecond_ && useSecantHessVec_ ) {
754  hist << " with " << ESecantToString(esec_) << " Hessian Approximation\n";
755  }
756  else {
757  hist << " with " << ESecantToString(esec_) << " Preconditioning and Hessian Approximation\n";
758  }
759  }
760  else {
761  hist << "\n";
762  }
763  if ( bndActive_ ) {
764  hist << "Trust-Region Model: " << ETrustRegionModelToString(TRmodel_) << "\n";
765  }
766  return hist.str();
767  }
768 
776  std::string print( AlgorithmState<Real> & algo_state, bool print_header = false ) const {
777  const Ptr<const StepState<Real>>& step_state = Step<Real>::getStepState();
778 
779  std::stringstream hist;
780  hist << std::scientific << std::setprecision(6);
781  if ( algo_state.iter == 0 ) {
782  hist << printName();
783  }
784  if ( print_header ) {
785  hist << printHeader();
786  }
787  if ( algo_state.iter == 0 ) {
788  hist << " ";
789  hist << std::setw(6) << std::left << algo_state.iter;
790  hist << std::setw(15) << std::left << algo_state.value;
791  hist << std::setw(15) << std::left << algo_state.gnorm;
792  hist << std::setw(15) << std::left << " ";
793  hist << std::setw(15) << std::left << step_state->searchSize;
794  hist << "\n";
795  }
796  else {
797  hist << " ";
798  hist << std::setw(6) << std::left << algo_state.iter;
799  hist << std::setw(15) << std::left << algo_state.value;
800  hist << std::setw(15) << std::left << algo_state.gnorm;
801  hist << std::setw(15) << std::left << algo_state.snorm;
802  hist << std::setw(15) << std::left << step_state->searchSize;
803  hist << std::setw(10) << std::left << algo_state.nfval;
804  hist << std::setw(10) << std::left << algo_state.ngrad;
805  hist << std::setw(10) << std::left << TRflag_;
807  hist << std::setw(10) << std::left << SPiter_;
808  hist << std::setw(10) << std::left << SPflag_;
809  }
810  hist << "\n";
811  }
812  return hist.str();
813  }
814 
815 }; // class Step
816 
817 } // namespace ROL
818 
819 #endif
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Contains definitions of custom data types in ROL.
Provides the interface to apply upper and lower bound constraints.
bool isActivated(void) const
Check if bounds are on.
virtual void projectInterior(Vector< Real > &x)
Project optimization variables into the interior of the feasible set.
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
void computeProjectedGradient(Vector< Real > &g, const Vector< Real > &x)
Compute projected gradient.
Provides the interface to evaluate objective functions.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:70
Provides the interface to compute optimization steps.
Definition: ROL_Step.hpp:68
ROL::Ptr< StepState< Real > > getState(void)
Definition: ROL_Step.hpp:73
const ROL::Ptr< const StepState< Real > > getStepState(void) const
Get state for step object.
Definition: ROL_Step.hpp:211
Provides the interface to compute optimization steps with trust regions.
TrustRegionStep(ROL::ParameterList &parlist)
Constructor.
int verbosity_
Print additional information to screen if > 0.
Ptr< Vector< Real > > gp_
Container for previous gradient vector.
Real beta_
Post-Smoothing rate for projected methods.
Real mu_
Post-Smoothing tolerance for projected methods.
Real scale0_
Scale for inexact gradient computation.
int SPiter_
Subproblem solver iteration count.
std::string print(AlgorithmState< Real > &algo_state, bool print_header=false) const
Print iterate status.
Real delMax_
Maximum trust-region radius.
void updateGradient(Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Update gradient to iteratively satisfy inexactness condition.
ETrustRegionModel TRmodel_
Trust-region subproblem model type.
Ptr< TrustRegion< Real > > trustRegion_
Container for trust-region solver object.
ESecant esec_
Secant type.
Ptr< TrustRegionModel< Real > > model_
Container for trust-region model.
TrustRegionStep(ROL::Ptr< Secant< Real > > &secant, ROL::ParameterList &parlist)
Constructor.
void initialize(Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Initialize step.
Ptr< Vector< Real > > xold_
Container for previous iteration vector.
bool useProjectedGrad_
Flag whether to use the projected gradient criticality measure.
bool bndActive_
Flag whether bound is activated.
ETrustRegion etr_
Trust-region subproblem solver type.
Ptr< Secant< Real > > secant_
Container for secant approximation.
std::string printHeader(void) const
Print iterate header.
int max_fval_
Maximum function evaluations in line-search for projected methods.
int SPflag_
Subproblem solver termination flag.
bool useSecantPrecond_
Flag whether to use a secant preconditioner.
Real computeCriticalityMeasure(const Vector< Real > &g, const Vector< Real > &x, BoundConstraint< Real > &bnd)
Compute the criticality measure.
ETrustRegionFlag TRflag_
Trust-region exit flag.
Real alpha_init_
Initial line-search parameter for projected methods.
Real scaleEps_
Scaling for epsilon-active sets.
void update(Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Update step, if successful.
std::vector< bool > useInexact_
Flags for inexact (0) objective function, (1) gradient, (2) Hessian.
bool useSecantHessVec_
Flag whether to use a secant Hessian.
void compute(Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Compute step.
Ptr< Vector< Real > > xnew_
Container for updated iteration vector.
std::string printName(void) const
Print step name.
Real scale1_
Scale for inexact gradient computation.
void parseParameterList(ROL::ParameterList &parlist)
Parse input ParameterList.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:84
virtual Real norm() const =0
Returns where .
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
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
ETrustRegionModel StringToETrustRegionModel(std::string s)
std::string NumberToString(T Number)
Definition: ROL_Types.hpp:81
ESecant StringToESecant(std::string s)
Definition: ROL_Types.hpp:541
@ SECANT_LBFGS
Definition: ROL_Types.hpp:485
ETrustRegion StringToETrustRegion(std::string s)
@ CG_FLAG_UNDEFINED
Definition: ROL_Types.hpp:825
@ CG_FLAG_SUCCESS
Definition: ROL_Types.hpp:820
bool isValidTrustRegionSubproblem(ETrustRegion etr, ETrustRegionModel etrm, bool isBnd)
std::string ETrustRegionModelToString(ETrustRegionModel tr)
@ TRUSTREGION_FLAG_POSPREDNEG
@ TRUSTREGION_FLAG_UNDEFINED
@ TRUSTREGION_FLAG_SUCCESS
std::string ESecantToString(ESecant tr)
Definition: ROL_Types.hpp:493
std::string ETrustRegionFlagToString(ETrustRegionFlag trf)
std::string ETrustRegionToString(ETrustRegion tr)
@ TRUSTREGION_MODEL_COLEMANLI
@ TRUSTREGION_MODEL_KELLEYSACHS
@ TRUSTREGION_MODEL_LINMORE
std::string ECGFlagToString(ECGFlag cgf)
Definition: ROL_Types.hpp:829
@ TRUSTREGION_DOUBLEDOGLEG
@ TRUSTREGION_TRUNCATEDCG
State for algorithm class. Will be used for restarts.
Definition: ROL_Types.hpp:143
ROL::Ptr< Vector< Real > > iterateVec
Definition: ROL_Types.hpp:157