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  obj.update(x,true,algo_state.iter);
494  }
495  // Build trust-region model
496  if (bnd.isActivated()) {
498  model_ = makePtr<KelleySachsModel<Real>>(obj,
499  bnd,
500  x,
501  *(step_state->gradientVec),
502  secant_,
505  }
506  else if ( TRmodel_ == TRUSTREGION_MODEL_COLEMANLI ) {
507  model_ = makePtr<ColemanLiModel<Real>>(obj,
508  bnd,
509  x,
510  *(step_state->gradientVec),
511  stepBackMax_,
514  secant_,
517  }
518  else if ( TRmodel_ == TRUSTREGION_MODEL_LINMORE ) {
519  model_ = makePtr<LinMoreModel<Real>>(obj,
520  bnd,
521  x,
522  *(step_state->gradientVec),
523  secant_,
526  }
527  else {
528  ROL_TEST_FOR_EXCEPTION( true, std::invalid_argument,
529  ">>> ERROR (TrustRegionStep): Invalid trust-region model!");
530  }
531  }
532  else {
533  model_ = makePtr<TrustRegionModel<Real>>(obj,
534  bnd,
535  x,
536  *(step_state->gradientVec),
537  secant_,
540  }
541  }
542 
553  void compute( Vector<Real> &s, const Vector<Real> &x,
555  AlgorithmState<Real> &algo_state ) {
556  // Get step state
557  Ptr<StepState<Real>> step_state = Step<Real>::getState();
558  // Build trust-region model
559  model_->update(obj,bnd,x,*step_state->gradientVec,secant_);
560  if (bnd.isActivated()) {
562 // Real eps = scaleEps_*algo_state.gnorm;
563  Real eps = scaleEps_ * std::min(std::pow(algo_state.gnorm,static_cast<Real>(0.75)),
564  static_cast<Real>(0.001));
565  dynamicPtrCast<KelleySachsModel<Real>>(model_)->setEpsilon(eps);
566  }
567  else if ( TRmodel_ == TRUSTREGION_MODEL_COLEMANLI ) {
568  dynamicPtrCast<ColemanLiModel<Real>>(model_)->setRadius(step_state->searchSize);
569  }
570  }
571  // Minimize trust-region model over trust-region constraint
572  SPflag_ = 0; SPiter_ = 0;
573  trustRegion_->run(s,algo_state.snorm,SPflag_,SPiter_,step_state->searchSize,*model_);
574  }
575 
588  const Vector<Real> &s,
589  Objective<Real> &obj,
591  AlgorithmState<Real> &algo_state ) {
592  // Get step state
593  Ptr<StepState<Real>> state = Step<Real>::getState();
594  // Store previous step for constraint computations
595  if ( bnd.isActivated() ) {
596  xold_->set(x);
597  }
598  // Update trust-region information;
599  // Performs a hard update on the objective function
601  state->nfval = 0;
602  state->ngrad = 0;
603  Real fold = algo_state.value;
604  Real fnew(0);
605  algo_state.iter++;
606  trustRegion_->update(x,fnew,state->searchSize,state->nfval,state->ngrad,TRflag_,
607  s,algo_state.snorm,fold,*(state->gradientVec),algo_state.iter,
608  obj,bnd,*model_);
609  algo_state.nfval += state->nfval;
610  algo_state.ngrad += state->ngrad;
611  state->flag = static_cast<int>(TRflag_);
612  state->SPiter = SPiter_;
613  state->SPflag = SPflag_;
614  // If step is accepted ...
615  // Compute new gradient and update secant storage
618  // Store previous gradient for secant update
620  gp_->set(*(state->gradientVec));
621  }
622  // Update objective function and approximate model
623  updateGradient(x,obj,bnd,algo_state);
624  // Update secant information
626  if ( bnd.isActivated() ) { // Compute new constrained step
627  xnew_->set(x);
628  xnew_->axpy(-static_cast<Real>(1),*xold_);
629  secant_->updateStorage(x,*(state->gradientVec),*gp_,*xnew_,algo_state.snorm,algo_state.iter+1);
630  }
631  else {
632  secant_->updateStorage(x,*(state->gradientVec),*gp_,s,algo_state.snorm,algo_state.iter+1);
633  }
634  }
635  // Update algorithm state
636  (algo_state.iterateVec)->set(x);
637  }
638  else {
639  if ( useInexact_[1] ) {
640  // Update objective function and approximate model
641  updateGradient(x,obj,bnd,algo_state);
642  }
643  }
644  // Update algorithm state
645  algo_state.value = fnew;
646  }
647 
652  std::string printHeader( void ) const {
653  std::stringstream hist;
654 
655  if(verbosity_>0) {
656  hist << std::string(114,'-') << "\n";
657 
658  hist << "Trust-Region status output definitions\n\n";
659 
660  hist << " iter - Number of iterates (steps taken) \n";
661  hist << " value - Objective function value \n";
662  hist << " gnorm - Norm of the gradient\n";
663  hist << " snorm - Norm of the step (update to optimization vector)\n";
664  hist << " delta - Trust-Region radius\n";
665  hist << " #fval - Number of times the objective function was evaluated\n";
666  hist << " #grad - Number of times the gradient was computed\n";
667 
668 
669 
670  hist << "\n";
671  hist << " tr_flag - Trust-Region flag" << "\n";
672  for( int flag = TRUSTREGION_FLAG_SUCCESS; flag != TRUSTREGION_FLAG_UNDEFINED; ++flag ) {
673  hist << " " << NumberToString(flag) << " - "
674  << ETrustRegionFlagToString(static_cast<ETrustRegionFlag>(flag)) << "\n";
675 
676  }
677 
678  if( etr_ == TRUSTREGION_TRUNCATEDCG ) {
679  hist << "\n";
680  hist << " iterCG - Number of Truncated CG iterations\n\n";
681  hist << " flagGC - Trust-Region Truncated CG flag" << "\n";
682  for( int flag = CG_FLAG_SUCCESS; flag != CG_FLAG_UNDEFINED; ++flag ) {
683  hist << " " << NumberToString(flag) << " - "
684  << ECGFlagToString(static_cast<ECGFlag>(flag)) << "\n";
685  }
686  }
687 
688  hist << std::string(114,'-') << "\n";
689  }
690 
691  hist << " ";
692  hist << std::setw(6) << std::left << "iter";
693  hist << std::setw(15) << std::left << "value";
694  hist << std::setw(15) << std::left << "gnorm";
695  hist << std::setw(15) << std::left << "snorm";
696  hist << std::setw(15) << std::left << "delta";
697  hist << std::setw(10) << std::left << "#fval";
698  hist << std::setw(10) << std::left << "#grad";
699  hist << std::setw(10) << std::left << "tr_flag";
701  hist << std::setw(10) << std::left << "iterCG";
702  hist << std::setw(10) << std::left << "flagCG";
703  }
704  hist << "\n";
705  return hist.str();
706  }
707 
712  std::string printName( void ) const {
713  std::stringstream hist;
714  hist << "\n" << ETrustRegionToString(etr_) << " Trust-Region Solver";
717  hist << " with " << ESecantToString(esec_) << " Preconditioning\n";
718  }
719  else if ( !useSecantPrecond_ && useSecantHessVec_ ) {
720  hist << " with " << ESecantToString(esec_) << " Hessian Approximation\n";
721  }
722  else {
723  hist << " with " << ESecantToString(esec_) << " Preconditioning and Hessian Approximation\n";
724  }
725  }
726  else {
727  hist << "\n";
728  }
729  if ( bndActive_ ) {
730  hist << "Trust-Region Model: " << ETrustRegionModelToString(TRmodel_) << "\n";
731  }
732  return hist.str();
733  }
734 
742  std::string print( AlgorithmState<Real> & algo_state, bool print_header = false ) const {
743  const Ptr<const StepState<Real>>& step_state = Step<Real>::getStepState();
744 
745  std::stringstream hist;
746  hist << std::scientific << std::setprecision(6);
747  if ( algo_state.iter == 0 ) {
748  hist << printName();
749  }
750  if ( print_header ) {
751  hist << printHeader();
752  }
753  if ( algo_state.iter == 0 ) {
754  hist << " ";
755  hist << std::setw(6) << std::left << algo_state.iter;
756  hist << std::setw(15) << std::left << algo_state.value;
757  hist << std::setw(15) << std::left << algo_state.gnorm;
758  hist << std::setw(15) << std::left << " ";
759  hist << std::setw(15) << std::left << step_state->searchSize;
760  hist << "\n";
761  }
762  else {
763  hist << " ";
764  hist << std::setw(6) << std::left << algo_state.iter;
765  hist << std::setw(15) << std::left << algo_state.value;
766  hist << std::setw(15) << std::left << algo_state.gnorm;
767  hist << std::setw(15) << std::left << algo_state.snorm;
768  hist << std::setw(15) << std::left << step_state->searchSize;
769  hist << std::setw(10) << std::left << algo_state.nfval;
770  hist << std::setw(10) << std::left << algo_state.ngrad;
771  hist << std::setw(10) << std::left << TRflag_;
773  hist << std::setw(10) << std::left << SPiter_;
774  hist << std::setw(10) << std::left << SPflag_;
775  }
776  hist << "\n";
777  }
778  return hist.str();
779  }
780 
781 }; // class Step
782 
783 } // namespace ROL
784 
785 #endif
std::string ECGFlagToString(ECGFlag cgf)
Definition: ROL_Types.hpp:829
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:68
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:73
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 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:211
Ptr< Vector< Real > > gp_
Container for previous gradient vector.