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_*std::min(algo_state.gnorm,state->searchSize);
251  Real gtol0 = gtol1 + one;
252  while ( gtol0 > gtol1 ) {
253  obj.gradient(*(state->gradientVec),x,gtol1);
254  algo_state.gnorm = computeCriticalityMeasure(*(state->gradientVec),x,bnd);
255  gtol0 = gtol1;
256  gtol1 = scale0_*std::min(algo_state.gnorm,state->searchSize);
257  }
258  algo_state.ngrad++;
259  }
260  else {
261  Real gtol = std::sqrt(ROL_EPSILON<Real>());
262  obj.gradient(*(state->gradientVec),x,gtol);
263  algo_state.ngrad++;
264  algo_state.gnorm = computeCriticalityMeasure(*(state->gradientVec),x,bnd);
265  }
266  }
267 
277  if ( bnd.isActivated() ) {
278  if ( useProjectedGrad_ ) {
279  gp_->set(g);
280  bnd.computeProjectedGradient( *gp_, x );
281  return gp_->norm();
282  }
283  else {
284  Real one(1);
285  xnew_->set(x);
286  xnew_->axpy(-one,g.dual());
287  bnd.project(*xnew_);
288  xnew_->axpy(-one,x);
289  return xnew_->norm();
290  }
291  }
292  else {
293  return g.norm();
294  }
295  }
296 
297 public:
298 
300  using Step<Real>::compute;
301  using Step<Real>::update;
302 
303  virtual ~TrustRegionStep() {}
304 
312  TrustRegionStep( ROL::ParameterList & parlist )
313  : Step<Real>(),
314  xnew_(nullPtr), xold_(nullPtr), gp_(nullPtr),
315  trustRegion_(nullPtr), model_(nullPtr),
318  SPflag_(0), SPiter_(0), bndActive_(false),
319  secant_(nullPtr), esec_(SECANT_LBFGS),
320  useSecantHessVec_(false), useSecantPrecond_(false),
321  scaleEps_(1), useProjectedGrad_(false),
322  alpha_init_(1), max_fval_(20), mu_(0.9999), beta_(0.01),
323  stepBackMax_(0.9999), stepBackScale_(1), singleReflect_(true),
324  scale0_(1), scale1_(1),
325  verbosity_(0) {
326  // Parse input parameterlist
327  parseParameterList(parlist);
328  // Create secant object
329  ROL::ParameterList &glist = parlist.sublist("General");
330  esec_ = StringToESecant(glist.sublist("Secant").get("Type","Limited-Memory BFGS"));
331  useSecantPrecond_ = glist.sublist("Secant").get("Use as Preconditioner", false);
332  useSecantHessVec_ = glist.sublist("Secant").get("Use as Hessian", false);
333  secant_ = SecantFactory<Real>(parlist);
334  }
335 
345  TrustRegionStep( ROL::Ptr<Secant<Real> > &secant, ROL::ParameterList &parlist )
346  : Step<Real>(),
347  xnew_(nullPtr), xold_(nullPtr), gp_(nullPtr),
348  trustRegion_(nullPtr), model_(nullPtr),
351  SPflag_(0), SPiter_(0), bndActive_(false),
352  secant_(nullPtr), esec_(SECANT_LBFGS),
353  useSecantHessVec_(false), useSecantPrecond_(false),
354  scaleEps_(1), useProjectedGrad_(false),
355  alpha_init_(1), max_fval_(20), mu_(0.9999), beta_(0.01),
356  stepBackMax_(0.9999), stepBackScale_(1), singleReflect_(true),
357  scale0_(1), scale1_(1),
358  verbosity_(0) {
359  // Parse input parameterlist
360  parseParameterList(parlist);
361  // Create secant object
362  ROL::ParameterList &glist = parlist.sublist("General");
363  useSecantPrecond_ = glist.sublist("Secant").get("Use as Preconditioner", false);
364  useSecantHessVec_ = glist.sublist("Secant").get("Use as Hessian", false);
365  if ( ROL::is_nullPtr(secant_) ) {
366  ROL::ParameterList Slist;
367  Slist.sublist("General").sublist("Secant").set("Type","Limited-Memory BFGS");
368  Slist.sublist("General").sublist("Secant").set("Maximum Storage",10);
369  secant_ = SecantFactory<Real>(Slist);
370  }
371  }
372 
381  void initialize( Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g,
383  AlgorithmState<Real> &algo_state ) {
385  throw Exception::NotImplemented(">>> ROL::TrustRegionStep : Invalid Trust Region Solver and Model pair!");
386  }
387  Real p1(0.1), oe10(1.e10), zero(0), one(1), half(0.5), three(3), two(2), six(6);
388  Ptr<StepState<Real>> step_state = Step<Real>::getState();
389  bndActive_ = bnd.isActivated();
390 
391  trustRegion_->initialize(x,s,g);
392 
393  Real htol = std::sqrt(ROL_EPSILON<Real>());
394  Real ftol = p1*ROL_OVERFLOW<Real>();
395 
396  step_state->descentVec = s.clone();
397  step_state->gradientVec = g.clone();
398 
399  if ( bnd.isActivated() ) {
400  // Make initial guess feasible
402  bnd.projectInterior(x);
403  }
404  else {
405  bnd.project(x);
406  }
407  xnew_ = x.clone();
408  xold_ = x.clone();
409  }
410  gp_ = g.clone();
411 
412  // Update approximate gradient and approximate objective function.
413  obj.update(x,true,algo_state.iter);
414  algo_state.snorm = oe10;
415  algo_state.value = obj.value(x,ftol);
416  algo_state.nfval++;
417  updateGradient(x,obj,bnd,algo_state);
418 
419  // Try to apply inverse Hessian
420  if ( !useSecantHessVec_ &&
422  try {
423  Ptr<Vector<Real>> v = g.clone();
424  Ptr<Vector<Real>> hv = x.clone();
425  obj.invHessVec(*hv,*v,x,htol);
426  }
427  catch (std::exception &e) {
428  useSecantHessVec_ = true;
429  }
430  }
431 
432  // Evaluate Objective Function at Cauchy Point
433  bool autoRad = false;
434  if ( step_state->searchSize <= zero ) {
435  autoRad = true;
436  Ptr<Vector<Real>> Bg = g.clone();
437  if ( useSecantHessVec_ ) {
438  secant_->applyB(*Bg,(step_state->gradientVec)->dual());
439  }
440  else {
441  obj.hessVec(*Bg,(step_state->gradientVec)->dual(),x,htol);
442  }
443  Real gBg = Bg->dot(*(step_state->gradientVec));
444  Real alpha = one;
445  if ( gBg > ROL_EPSILON<Real>() ) {
446  alpha = algo_state.gnorm*algo_state.gnorm/gBg;
447  }
448  // Evaluate the objective function at the Cauchy point
449  Ptr<Vector<Real>> cp = s.clone();
450  cp->set((step_state->gradientVec)->dual());
451  cp->scale(-alpha);
452  Ptr<Vector<Real>> xcp = x.clone();
453  xcp->set(x);
454  xcp->plus(*cp);
455  if ( bnd.isActivated() ) {
456  bnd.project(*xcp);
457  }
458  obj.update(*xcp);
459  Real fnew = obj.value(*xcp,ftol); // MUST DO SOMETHING HERE WITH FTOL
460  algo_state.nfval++;
461  // Perform cubic interpolation to determine initial trust region radius
462  Real gs = cp->dot((step_state->gradientVec)->dual());
463  Real a = fnew - algo_state.value - gs - half*alpha*alpha*gBg;
464  if ( std::abs(a) < ROL_EPSILON<Real>() ) {
465  // a = 0 implies the objective is quadratic in the negative gradient direction
466  step_state->searchSize = std::min(alpha*algo_state.gnorm,delMax_);
467  }
468  else {
469  Real b = half*alpha*alpha*gBg;
470  Real c = gs;
471  if ( b*b-three*a*c > ROL_EPSILON<Real>() ) {
472  // There is at least one critical point
473  Real t1 = (-b-std::sqrt(b*b-three*a*c))/(three*a);
474  Real t2 = (-b+std::sqrt(b*b-three*a*c))/(three*a);
475  if ( six*a*t1 + two*b > zero ) {
476  // t1 is the minimizer
477  step_state->searchSize = std::min(t1*alpha*algo_state.gnorm,delMax_);
478  }
479  else {
480  // t2 is the minimizer
481  step_state->searchSize = std::min(t2*alpha*algo_state.gnorm,delMax_);
482  }
483  }
484  else {
485  step_state->searchSize = std::min(alpha*algo_state.gnorm,delMax_);
486  }
487  }
488  if (step_state->searchSize <= ROL_EPSILON<Real>()*algo_state.gnorm && autoRad) {
489  step_state->searchSize = one;
490  }
491  }
492  // Build trust-region model
493  if (bnd.isActivated()) {
495  model_ = makePtr<KelleySachsModel<Real>>(obj,
496  bnd,
497  x,
498  *(step_state->gradientVec),
499  secant_,
502  }
503  else if ( TRmodel_ == TRUSTREGION_MODEL_COLEMANLI ) {
504  model_ = makePtr<ColemanLiModel<Real>>(obj,
505  bnd,
506  x,
507  *(step_state->gradientVec),
508  stepBackMax_,
511  secant_,
514  }
515  else if ( TRmodel_ == TRUSTREGION_MODEL_LINMORE ) {
516  model_ = makePtr<LinMoreModel<Real>>(obj,
517  bnd,
518  x,
519  *(step_state->gradientVec),
520  secant_,
523  }
524  else {
525  ROL_TEST_FOR_EXCEPTION( true, std::invalid_argument,
526  ">>> ERROR (TrustRegionStep): Invalid trust-region model!");
527  }
528  }
529  else {
530  model_ = makePtr<TrustRegionModel<Real>>(obj,
531  bnd,
532  x,
533  *(step_state->gradientVec),
534  secant_,
537  }
538  }
539 
550  void compute( Vector<Real> &s, const Vector<Real> &x,
552  AlgorithmState<Real> &algo_state ) {
553  // Get step state
554  Ptr<StepState<Real>> step_state = Step<Real>::getState();
555  // Build trust-region model
556  model_->update(obj,bnd,x,*step_state->gradientVec,secant_);
557  if (bnd.isActivated()) {
559 // Real eps = scaleEps_*algo_state.gnorm;
560  Real eps = scaleEps_ * std::min(std::pow(algo_state.gnorm,static_cast<Real>(0.75)),
561  static_cast<Real>(0.001));
562  dynamicPtrCast<KelleySachsModel<Real>>(model_)->setEpsilon(eps);
563  }
564  else if ( TRmodel_ == TRUSTREGION_MODEL_COLEMANLI ) {
565  dynamicPtrCast<ColemanLiModel<Real>>(model_)->setRadius(step_state->searchSize);
566  }
567  }
568  // Minimize trust-region model over trust-region constraint
569  SPflag_ = 0; SPiter_ = 0;
570  trustRegion_->run(s,algo_state.snorm,SPflag_,SPiter_,step_state->searchSize,*model_);
571  }
572 
585  const Vector<Real> &s,
586  Objective<Real> &obj,
588  AlgorithmState<Real> &algo_state ) {
589  // Get step state
590  Ptr<StepState<Real>> state = Step<Real>::getState();
591  // Store previous step for constraint computations
592  if ( bnd.isActivated() ) {
593  xold_->set(x);
594  }
595  // Update trust-region information;
596  // Performs a hard update on the objective function
598  state->nfval = 0;
599  state->ngrad = 0;
600  Real fold = algo_state.value;
601  Real fnew(0);
602  algo_state.iter++;
603  trustRegion_->update(x,fnew,state->searchSize,state->nfval,state->ngrad,TRflag_,
604  s,algo_state.snorm,fold,*(state->gradientVec),algo_state.iter,
605  obj,bnd,*model_);
606  algo_state.nfval += state->nfval;
607  algo_state.ngrad += state->ngrad;
608  state->flag = static_cast<int>(TRflag_);
609  state->SPiter = SPiter_;
610  state->SPflag = SPflag_;
611  // If step is accepted ...
612  // Compute new gradient and update secant storage
615  // Perform line search (smoothing) to ensure decrease
617  Real tol = std::sqrt(ROL_EPSILON<Real>());
618  // Compute new gradient
619  obj.gradient(*gp_,x,tol); // MUST DO SOMETHING HERE WITH TOL
620  algo_state.ngrad++;
621  // Compute smoothed step
622  Real alpha(1);
623  xnew_->set(x);
624  xnew_->axpy(-alpha*alpha_init_,gp_->dual());
625  bnd.project(*xnew_);
626  // Compute new objective value
627  obj.update(*xnew_,true,algo_state.iter);
628  Real ftmp = obj.value(*xnew_,tol); // MUST DO SOMETHING HERE WITH TOL
629  algo_state.nfval++;
630  // Perform smoothing
631  int cnt = 0;
632  alpha = static_cast<Real>(1)/alpha_init_;
633  while ( (fnew-ftmp) <= mu_*(fnew-fold) ) {
634  xnew_->set(x);
635  xnew_->axpy(-alpha*alpha_init_,gp_->dual());
636  bnd.project(*xnew_);
637  obj.update(*xnew_,true,algo_state.iter);
638  ftmp = obj.value(*xnew_,tol); // MUST DO SOMETHING HERE WITH TOL
639  algo_state.nfval++;
640  if ( cnt >= max_fval_ ) {
641  break;
642  }
643  alpha *= beta_;
644  cnt++;
645  }
646  // Store objective function and iteration information
647  fnew = ftmp;
648  x.set(*xnew_);
649  }
650  // Store previous gradient for secant update
652  gp_->set(*(state->gradientVec));
653  }
654  // Update objective function and approximate model
655  updateGradient(x,obj,bnd,algo_state);
656  // Update secant information
658  if ( bnd.isActivated() ) { // Compute new constrained step
659  xnew_->set(x);
660  xnew_->axpy(-static_cast<Real>(1),*xold_);
661  secant_->updateStorage(x,*(state->gradientVec),*gp_,*xnew_,algo_state.snorm,algo_state.iter+1);
662  }
663  else {
664  secant_->updateStorage(x,*(state->gradientVec),*gp_,s,algo_state.snorm,algo_state.iter+1);
665  }
666  }
667  // Update algorithm state
668  (algo_state.iterateVec)->set(x);
669  }
670  else {
671  if ( useInexact_[1] ) {
672  // Update objective function and approximate model
673  updateGradient(x,obj,bnd,algo_state);
674  }
675  }
676  // Update algorithm state
677  algo_state.value = fnew;
678  }
679 
684  std::string printHeader( void ) const {
685  std::stringstream hist;
686 
687  if(verbosity_>0) {
688  hist << std::string(114,'-') << "\n";
689 
690  hist << "Trust-Region status output definitions\n\n";
691 
692  hist << " iter - Number of iterates (steps taken) \n";
693  hist << " value - Objective function value \n";
694  hist << " gnorm - Norm of the gradient\n";
695  hist << " snorm - Norm of the step (update to optimization vector)\n";
696  hist << " delta - Trust-Region radius\n";
697  hist << " #fval - Number of times the objective function was evaluated\n";
698  hist << " #grad - Number of times the gradient was computed\n";
699 
700 
701 
702  hist << "\n";
703  hist << " tr_flag - Trust-Region flag" << "\n";
704  for( int flag = TRUSTREGION_FLAG_SUCCESS; flag != TRUSTREGION_FLAG_UNDEFINED; ++flag ) {
705  hist << " " << NumberToString(flag) << " - "
706  << ETrustRegionFlagToString(static_cast<ETrustRegionFlag>(flag)) << "\n";
707 
708  }
709 
710  if( etr_ == TRUSTREGION_TRUNCATEDCG ) {
711  hist << "\n";
712  hist << " iterCG - Number of Truncated CG iterations\n\n";
713  hist << " flagGC - Trust-Region Truncated CG flag" << "\n";
714  for( int flag = CG_FLAG_SUCCESS; flag != CG_FLAG_UNDEFINED; ++flag ) {
715  hist << " " << NumberToString(flag) << " - "
716  << ECGFlagToString(static_cast<ECGFlag>(flag)) << "\n";
717  }
718  }
719 
720  hist << std::string(114,'-') << "\n";
721  }
722 
723  hist << " ";
724  hist << std::setw(6) << std::left << "iter";
725  hist << std::setw(15) << std::left << "value";
726  hist << std::setw(15) << std::left << "gnorm";
727  hist << std::setw(15) << std::left << "snorm";
728  hist << std::setw(15) << std::left << "delta";
729  hist << std::setw(10) << std::left << "#fval";
730  hist << std::setw(10) << std::left << "#grad";
731  hist << std::setw(10) << std::left << "tr_flag";
733  hist << std::setw(10) << std::left << "iterCG";
734  hist << std::setw(10) << std::left << "flagCG";
735  }
736  hist << "\n";
737  return hist.str();
738  }
739 
744  std::string printName( void ) const {
745  std::stringstream hist;
746  hist << "\n" << ETrustRegionToString(etr_) << " Trust-Region Solver";
749  hist << " with " << ESecantToString(esec_) << " Preconditioning\n";
750  }
751  else if ( !useSecantPrecond_ && useSecantHessVec_ ) {
752  hist << " with " << ESecantToString(esec_) << " Hessian Approximation\n";
753  }
754  else {
755  hist << " with " << ESecantToString(esec_) << " Preconditioning and Hessian Approximation\n";
756  }
757  }
758  else {
759  hist << "\n";
760  }
761  if ( bndActive_ ) {
762  hist << "Trust-Region Model: " << ETrustRegionModelToString(TRmodel_) << "\n";
763  }
764  return hist.str();
765  }
766 
774  std::string print( AlgorithmState<Real> & algo_state, bool print_header = false ) const {
775  const Ptr<const StepState<Real>>& step_state = Step<Real>::getStepState();
776 
777  std::stringstream hist;
778  hist << std::scientific << std::setprecision(6);
779  if ( algo_state.iter == 0 ) {
780  hist << printName();
781  }
782  if ( print_header ) {
783  hist << printHeader();
784  }
785  if ( algo_state.iter == 0 ) {
786  hist << " ";
787  hist << std::setw(6) << std::left << algo_state.iter;
788  hist << std::setw(15) << std::left << algo_state.value;
789  hist << std::setw(15) << std::left << algo_state.gnorm;
790  hist << std::setw(15) << std::left << " ";
791  hist << std::setw(15) << std::left << step_state->searchSize;
792  hist << "\n";
793  }
794  else {
795  hist << " ";
796  hist << std::setw(6) << std::left << algo_state.iter;
797  hist << std::setw(15) << std::left << algo_state.value;
798  hist << std::setw(15) << std::left << algo_state.gnorm;
799  hist << std::setw(15) << std::left << algo_state.snorm;
800  hist << std::setw(15) << std::left << step_state->searchSize;
801  hist << std::setw(10) << std::left << algo_state.nfval;
802  hist << std::setw(10) << std::left << algo_state.ngrad;
803  hist << std::setw(10) << std::left << TRflag_;
805  hist << std::setw(10) << std::left << SPiter_;
806  hist << std::setw(10) << std::left << SPflag_;
807  }
808  hist << "\n";
809  }
810  return hist.str();
811  }
812 
813 }; // class Step
814 
815 } // namespace ROL
816 
817 #endif
std::string ECGFlagToString(ECGFlag cgf)
Definition: ROL_Types.hpp:900
Provides the interface to evaluate objective functions.
ESecant esec_
Secant type.
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
bool bndActive_
Flag whether bound is activated.
std::string ETrustRegionModelToString(ETrustRegionModel tr)
Real mu_
Post-Smoothing tolerance for projected methods.
virtual void projectInterior(Vector< Real > &x)
Project optimization variables into the interior of the feasible set.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
bool useSecantPrecond_
Flag whether to use a secant preconditioner.
bool isActivated(void) const
Check if bounds are on.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
Provides the interface to compute optimization steps.
Definition: ROL_Step.hpp:69
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.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Contains definitions of custom data types in ROL.
bool isValidTrustRegionSubproblem(ETrustRegion etr, ETrustRegionModel etrm, bool isBnd)
Real alpha_init_
Initial line-search parameter for projected methods.
ESecant StringToESecant(std::string s)
Definition: ROL_Types.hpp:541
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
std::string printName(void) const
Print step name.
TrustRegionStep(ROL::ParameterList &parlist)
Constructor.
Real computeCriticalityMeasure(const Vector< Real > &g, const Vector< Real > &x, BoundConstraint< Real > &bnd)
Compute the criticality measure.
ETrustRegionModel TRmodel_
Trust-region subproblem model type.
State for algorithm class. Will be used for restarts.
Definition: ROL_Types.hpp:143
Real scaleEps_
Scaling for epsilon-active sets.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
std::string NumberToString(T Number)
Definition: ROL_Types.hpp:81
ETrustRegion etr_
Trust-region subproblem solver type.
std::string printHeader(void) const
Print iterate header.
void compute(Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Compute step.
ESecant
Enumeration of secant update algorithms.
Definition: ROL_Types.hpp:484
ETrustRegionModel StringToETrustRegionModel(std::string s)
int SPflag_
Subproblem solver termination flag.
ROL::Ptr< StepState< Real > > getState(void)
Definition: ROL_Step.hpp:74
ETrustRegionFlag TRflag_
Trust-region exit flag.
bool useSecantHessVec_
Flag whether to use a secant Hessian.
int SPiter_
Subproblem solver iteration count.
std::vector< bool > useInexact_
Flags for inexact (0) objective function, (1) gradient, (2) Hessian.
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:70
Real scale0_
Scale for inexact gradient computation.
Ptr< TrustRegionModel< Real > > model_
Container for trust-region model.
ROL::Ptr< Vector< Real > > iterateVec
Definition: ROL_Types.hpp:157
Real delMax_
Maximum trust-region radius.
virtual void invHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply inverse Hessian approximation to vector.
Ptr< Vector< Real > > xnew_
Container for updated iteration vector.
std::string print(AlgorithmState< Real > &algo_state, bool print_header=false) const
Print iterate status.
ETrustRegionModel
Enumeration of trust-region model types.
Provides the interface to apply upper and lower bound constraints.
void computeProjectedGradient(Vector< Real > &g, const Vector< Real > &x)
Compute projected gradient.
int verbosity_
Print additional information to screen if &gt; 0.
Ptr< TrustRegion< Real > > trustRegion_
Container for trust-region solver object.
ETrustRegion StringToETrustRegion(std::string s)
Ptr< Vector< Real > > xold_
Container for previous iteration vector.
void update(Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Update step, if successful.
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
virtual Real norm() const =0
Returns where .
int max_fval_
Maximum function evaluations in line-search for projected methods.
Ptr< Secant< Real > > secant_
Container for secant approximation.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
ETrustRegion
Enumeration of trust-region solver types.
std::string ETrustRegionFlagToString(ETrustRegionFlag trf)
ETrustRegionFlag
Enumation of flags used by trust-region solvers.
std::string ETrustRegionToString(ETrustRegion tr)
std::string ESecantToString(ESecant tr)
Definition: ROL_Types.hpp:493
bool useProjectedGrad_
Flag whether to use the projected gradient criticality measure.
Real beta_
Post-Smoothing rate for projected methods.
Real scale1_
Scale for inexact gradient computation.
void parseParameterList(ROL::ParameterList &parlist)
Parse input ParameterList.
void updateGradient(Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Update gradient to iteratively satisfy inexactness condition.
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
Provides the interface to compute optimization steps with trust regions.
const ROL::Ptr< const StepState< Real > > getStepState(void) const
Get state for step object.
Definition: ROL_Step.hpp:293
Ptr< Vector< Real > > gp_
Container for previous gradient vector.