44 #ifndef ROL_BUNDLE_STEP_H 45 #define ROL_BUNDLE_STEP_H 56 #include "Teuchos_ParameterList.hpp" 57 #include "Teuchos_RCP.hpp" 82 Teuchos::RCP<Vector<Real> >
y_;
118 : bundle_(Teuchos::null), lineSearch_(Teuchos::null),
119 QPiter_(0), QPmaxit_(0), QPtol_(0), step_flag_(0),
120 y_(Teuchos::null), linErrNew_(0), valueNew_(0),
121 aggSubGradNew_(Teuchos::null), aggSubGradOldNorm_(0),
122 aggLinErrNew_(0), aggLinErrOld_(0), aggDistMeasNew_(0),
123 T_(0), tol_(0), m1_(0), m2_(0), m3_(0), nu_(0),
124 ls_maxit_(0), first_print_(true), isConvex_(false),
126 Real zero(0), 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);
128 state->searchSize = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Initial Trust-Region Parameter", oe3);
129 T_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Maximum Trust-Region Parameter", oe8);
130 tol_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Epsilon Solution Tolerance", oem6);
131 m1_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Upper Threshold for Serious Step", p1);
132 m2_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Lower Threshold for Serious Step", p2);
133 m3_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Upper Threshold for Null Step", p9);
134 nu_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Tolerance for Trust-Region Parameter", oem3);
137 Real coeff = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Distance Measure Coefficient", zero);
138 unsigned maxSize = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Maximum Bundle Size", 200);
139 unsigned remSize = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Removal Size for Bundle Update", 2);
140 if ( parlist.sublist(
"Step").sublist(
"Bundle").get(
"Cutting Plane Solver",0) == 1 ) {
142 bundle_ = Teuchos::rcp(
new Bundle<Real>(maxSize,coeff,remSize));
145 bundle_ = Teuchos::rcp(
new Bundle<Real>(maxSize,coeff,remSize));
147 isConvex_ = ((coeff == zero) ?
true :
false);
150 QPtol_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Cutting Plane Tolerance", oem8);
151 QPmaxit_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Cutting Plane Iteration Limit", 1000);
155 = parlist.sublist(
"Step").sublist(
"Line Search").get(
"Maximum Number of Function Evaluations",20);
157 lineSearch_ = LineSearchFactory<Real>(parlist);
166 Real searchSize = state->searchSize;
168 state->searchSize = searchSize;
170 bundle_->initialize(*(state->gradientVec));
174 aggSubGradNew_ = g.
clone();
175 aggSubGradOldNorm_ = algo_state.
gnorm;
178 lineSearch_->initialize(x,x,g,obj,con);
185 first_print_ =
false;
186 QPiter_ = (step_flag_ ? 0 :
QPiter_);
187 Real v(0), l(0), u =
T_, gd(0);
188 Real zero(0), two(2), half(0.5);
194 QPiter_ += bundle_->solveDual(state->searchSize,QPmaxit_,QPtol_);
195 bundle_->aggregate(*aggSubGradNew_,aggLinErrNew_,aggDistMeasNew_);
201 s.
set(aggSubGradNew_->dual()); s.
scale(-state->searchSize);
211 algo_state.
flag =
true;
216 y_->set(x); y_->plus(s);
218 valueNew_ = obj.
value(*y_,ftol_);
220 obj.
gradient(*(state->gradientVec),*y_,ftol_);
223 gd = s.
dot(state->gradientVec->dual());
224 linErrNew_ = algo_state.
value - (valueNew_ - gd);
226 bool SS1 = (valueNew_-algo_state.
value < m1_*v);
228 bool NS2a = (bundle_->computeAlpha(algo_state.
snorm,linErrNew_) <= m3_*
aggLinErrOld_);
229 bool NS2b = (std::abs(algo_state.
value-valueNew_) <= aggSubGradOldNorm_ +
aggLinErrOld_);
233 bool SS2 = (gd >= m2_*v || state->searchSize >= T_-
nu_);
240 l = state->searchSize;
241 state->searchSize = half*(u+l);
245 if ( NS2a || NS2b ) {
252 u = state->searchSize;
253 state->searchSize = half*(u+l);
259 bool NS3 = (gd - bundle_->computeAlpha(algo_state.
snorm,linErrNew_) >= m2_*v);
266 if ( NS2a || NS2b ) {
276 int ls_nfval = 0, ls_ngrad = 0;
277 lineSearch_->run(alpha,valueNew_,ls_nfval,ls_ngrad,gd,s,x,obj,con);
278 if ( ls_nfval == ls_maxit_ ) {
292 u = state->searchSize;
293 state->searchSize = half*(u+l);
298 u = state->searchSize;
299 state->searchSize = half*(u+l);
316 if ( !algo_state.
flag ) {
320 bundle_->reset(*aggSubGradNew_,aggLinErrNew_,algo_state.
snorm);
327 Real valueOld = algo_state.
value;
329 bundle_->update(step_flag_,valueNew_-valueOld,algo_state.
snorm,*(state->gradientVec),s);
333 bundle_->update(step_flag_,linErrNew_,algo_state.
snorm,*(state->gradientVec),s);
340 algo_state.
gnorm = (state->gradientVec)->norm();
347 std::stringstream hist;
349 hist << std::setw(6) << std::left <<
"iter";
350 hist << std::setw(15) << std::left <<
"value";
351 hist << std::setw(15) << std::left <<
"gnorm";
352 hist << std::setw(15) << std::left <<
"snorm";
353 hist << std::setw(10) << std::left <<
"#fval";
354 hist << std::setw(10) << std::left <<
"#grad";
355 hist << std::setw(15) << std::left <<
"znorm";
356 hist << std::setw(15) << std::left <<
"alpha";
357 hist << std::setw(15) << std::left <<
"TRparam";
358 hist << std::setw(10) << std::left <<
"QPiter";
364 std::stringstream hist;
365 hist <<
"\n" <<
"Bundle Trust-Region Algorithm \n";
371 std::stringstream hist;
372 hist << std::scientific << std::setprecision(6);
373 if ( algo_state.
iter == 0 && first_print_ ) {
375 if ( print_header ) {
379 hist << std::setw(6) << std::left << algo_state.
iter;
380 hist << std::setw(15) << std::left << algo_state.
value;
381 hist << std::setw(15) << std::left << algo_state.
gnorm;
384 if ( step_flag_ && algo_state.
iter > 0 ) {
385 if ( print_header ) {
390 hist << std::setw(6) << std::left << algo_state.
iter;
391 hist << std::setw(15) << std::left << algo_state.
value;
392 hist << std::setw(15) << std::left << algo_state.
gnorm;
393 hist << std::setw(15) << std::left << algo_state.
snorm;
394 hist << std::setw(10) << std::left << algo_state.
nfval;
395 hist << std::setw(10) << std::left << algo_state.
ngrad;
398 hist << std::setw(15) << std::left << state->searchSize;
399 hist << std::setw(10) << std::left <<
QPiter_;
Provides the interface to evaluate objective functions.
virtual void scale(const Real alpha)=0
Compute where .
virtual void plus(const Vector &x)=0
Compute , where .
std::string printHeader(void) const
Print iterate header.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
Provides the interface to compute optimization steps.
Real aggregateGradientNorm
Teuchos::RCP< StepState< Real > > getState(void)
Contains definitions of custom data types in ROL.
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
virtual Real dot(const Vector &x) const =0
Compute where .
State for algorithm class. Will be used for restarts.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
const Teuchos::RCP< const StepState< Real > > getStepState(void) const
Get state for step object.
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Initialize step with bound constraint.
Teuchos::RCP< Vector< Real > > y_
std::string print(AlgorithmState< Real > &algo_state, bool print_header=false) const
Print iterate status.
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.
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.
Teuchos::RCP< Bundle< Real > > bundle_
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.
Teuchos::RCP< Vector< Real > > iterateVec
virtual void set(const Vector &x)
Set where .
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
Teuchos::RCP< Vector< Real > > aggSubGradNew_
std::string printName(void) const
Print step name.
Teuchos::RCP< LineSearch< Real > > lineSearch_
BundleStep(Teuchos::ParameterList &parlist)
Provides the interface for and implments a bundle.