10 #ifndef ROL_BUNDLE_STEP_H
11 #define ROL_BUNDLE_STEP_H
22 #include "ROL_ParameterList.hpp"
23 #include "ROL_Ptr.hpp"
48 ROL::Ptr<Vector<Real> >
y_;
94 Real
zero(0), two(2), oem3(1.e-3), oem6(1.e-6), oem8(1.e-8), p1(0.1), p2(0.2), p9(0.9), oe3(1.e3), oe8(1.e8);
96 state->searchSize = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Initial Trust-Region Parameter", oe3);
97 T_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Maximum Trust-Region Parameter", oe8);
98 tol_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Epsilon Solution Tolerance", oem6);
99 m1_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Upper Threshold for Serious Step", p1);
100 m2_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Lower Threshold for Serious Step", p2);
101 m3_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Upper Threshold for Null Step", p9);
102 nu_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Tolerance for Trust-Region Parameter", oem3);
105 Real coeff = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Distance Measure Coefficient",
zero);
106 Real omega = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Locality Measure Coefficient", two);
107 unsigned maxSize = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Maximum Bundle Size", 200);
108 unsigned remSize = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Removal Size for Bundle Update", 2);
109 int cps = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Cutting Plane Solver",0);
111 bundle_ = ROL::makePtr<Bundle_TT<Real>>(maxSize,coeff,omega,remSize);
115 bundle_ = ROL::makePtr<Bundle_AS<Real>>(maxSize,coeff,omega,remSize);
120 QPtol_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Cutting Plane Tolerance", oem8);
121 QPmaxit_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Cutting Plane Iteration Limit", 1000);
125 = parlist.sublist(
"Step").sublist(
"Line Search").get(
"Maximum Number of Function Evaluations",20);
131 verbosity_ = parlist.sublist(
"General").get(
"Print Verbosity", 0);
139 Real searchSize = state->searchSize;
141 state->searchSize = searchSize;
143 bundle_->initialize(*(state->gradientVec));
160 Real v(0), l(0), u =
T_, gd(0);
161 Real
zero(0), two(2), half(0.5);
171 std::cout << std::endl;
172 std::cout <<
" Computation of aggregrate quantities" << std::endl;
174 std::cout <<
" Aggregate linearization error: " <<
aggLinErrNew_ << std::endl;
175 std::cout <<
" Aggregate distance measure: " <<
aggDistMeasNew_ << std::endl;
184 std::cout << std::endl;
185 std::cout <<
" Solve cutting plan subproblem" << std::endl;
186 std::cout <<
" Cutting plan objective value: " << v << std::endl;
187 std::cout <<
" Norm of computed step: " << algo_state.
snorm << std::endl;
188 std::cout <<
" 'Trust-region' radius: " << state->searchSize << std::endl;
198 algo_state.
flag =
true;
207 algo_state.
flag =
true;
211 y_->set(x);
y_->plus(s);
218 gd = s.
dot(state->gradientVec->dual());
221 Real eps =
static_cast<Real
>(10)*ROL_EPSILON<Real>();
222 Real del = eps*std::max(static_cast<Real>(1),std::abs(algo_state.
value));
226 if (std::abs(Df) < eps && std::abs(Dm) < eps) {
237 std::cout << std::endl;
238 std::cout <<
" Check for serious/null step" << std::endl;
239 std::cout <<
" Serious step test SS(i): " << SS1 << std::endl;
240 std::cout <<
" -> Left hand side: " <<
valueNew_-algo_state.
value << std::endl;
241 std::cout <<
" -> Right hand side: " <<
m1_*v << std::endl;
242 std::cout <<
" Null step test NS(iia): " << NS2a << std::endl;
245 std::cout <<
" Null step test NS(iib): " << NS2b << std::endl;
246 std::cout <<
" -> Left hand side: " << std::abs(algo_state.
value-
valueNew_) << std::endl;
252 bool SS2 = (gd >=
m2_*v || state->searchSize >=
T_-
nu_);
254 std::cout <<
" Serious step test SS(iia): " << (gd >=
m2_*v) << std::endl;
255 std::cout <<
" -> Left hand side: " << gd << std::endl;
256 std::cout <<
" -> Right hand side: " <<
m2_*v << std::endl;
257 std::cout <<
" Serious step test SS(iia): " << (state->searchSize >=
T_-
nu_) << std::endl;
258 std::cout <<
" -> Left hand side: " << state->searchSize << std::endl;
259 std::cout <<
" -> Right hand side: " <<
T_-
nu_ << std::endl;
265 std::cout <<
" Serious step taken" << std::endl;
270 l = state->searchSize;
271 state->searchSize = half*(u+l);
273 std::cout <<
" Increase 'trust-region' radius: " << state->searchSize << std::endl;
278 if ( NS2a || NS2b ) {
283 std::cout <<
" Null step taken" << std::endl;
288 u = state->searchSize;
289 state->searchSize = half*(u+l);
291 std::cout <<
" Decrease 'trust-region' radius: " << state->searchSize << std::endl;
300 std::cout <<
" Null step test NS(iii): " << NS3 << std::endl;
301 std::cout <<
" -> Left hand side: " << gd -
bundle_->computeAlpha(algo_state.
snorm,
linErrNew_) << std::endl;
302 std::cout <<
" -> Right hand side: " <<
m2_*v << std::endl;
310 if ( NS2a || NS2b ) {
320 int ls_nfval = 0, ls_ngrad = 0;
336 u = state->searchSize;
337 state->searchSize = half*(u+l);
342 u = state->searchSize;
343 state->searchSize = half*(u+l);
362 if ( !algo_state.
flag ) {
373 Real valueOld = algo_state.
value;
386 algo_state.
gnorm = (state->gradientVec)->norm();
393 std::stringstream hist;
395 hist << std::setw(6) << std::left <<
"iter";
396 hist << std::setw(15) << std::left <<
"value";
397 hist << std::setw(15) << std::left <<
"gnorm";
398 hist << std::setw(15) << std::left <<
"snorm";
399 hist << std::setw(10) << std::left <<
"#fval";
400 hist << std::setw(10) << std::left <<
"#grad";
401 hist << std::setw(15) << std::left <<
"znorm";
402 hist << std::setw(15) << std::left <<
"alpha";
403 hist << std::setw(15) << std::left <<
"TRparam";
404 hist << std::setw(10) << std::left <<
"QPiter";
410 std::stringstream hist;
411 hist <<
"\n" <<
"Bundle Trust-Region Algorithm \n";
417 std::stringstream hist;
418 hist << std::scientific << std::setprecision(6);
421 if ( print_header ) {
425 hist << std::setw(6) << std::left << algo_state.
iter;
426 hist << std::setw(15) << std::left << algo_state.
value;
427 hist << std::setw(15) << std::left << algo_state.
gnorm;
431 if ( print_header ) {
436 hist << std::setw(6) << std::left << algo_state.
iter;
437 hist << std::setw(15) << std::left << algo_state.
value;
438 hist << std::setw(15) << std::left << algo_state.
gnorm;
439 hist << std::setw(15) << std::left << algo_state.
snorm;
440 hist << std::setw(10) << std::left << algo_state.
nfval;
441 hist << std::setw(10) << std::left << algo_state.
ngrad;
444 hist << std::setw(15) << std::left << state->searchSize;
445 hist << std::setw(10) << std::left <<
QPiter_;
Provides the interface to evaluate objective functions.
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
ROL::Ptr< Vector< Real > > y_
virtual void plus(const Vector &x)=0
Compute , where .
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
Provides the interface to compute optimization steps.
Real aggregateGradientNorm
Contains definitions of custom data types in ROL.
BundleStep(ROL::ParameterList &parlist)
std::string printName(void) const
Print step name.
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
virtual Real dot(const Vector &x) const =0
Compute where .
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
State for algorithm class. Will be used for restarts.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
ROL::Ptr< StepState< Real > > getState(void)
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Initialize step with bound constraint.
ROL::Ptr< Vector< Real > > iterateVec
std::string printHeader(void) const
Print iterate header.
Provides the interface to apply upper and lower bound constraints.
void update(Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Update step, if successful.
std::string print(AlgorithmState< Real > &algo_state, bool print_header=false) const
Print iterate status.
void compute(Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Compute step.
Provides the interface to compute bundle trust-region steps.
ROL::Ptr< Vector< Real > > aggSubGradNew_
virtual void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Initialize step with bound constraint.
virtual void set(const Vector &x)
Set where .
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
ROL::Ptr< Bundle< Real > > bundle_
ROL::Ptr< LineSearch< Real > > lineSearch_
const ROL::Ptr< const StepState< Real > > getStepState(void) const
Get state for step object.