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_;
120 : bundle_(Teuchos::null), lineSearch_(Teuchos::null),
121 QPiter_(0), QPmaxit_(0), QPtol_(0), step_flag_(0),
122 y_(Teuchos::null), linErrNew_(0), valueNew_(0),
123 aggSubGradNew_(Teuchos::null), aggSubGradOldNorm_(0),
124 aggLinErrNew_(0), aggLinErrOld_(0), aggDistMeasNew_(0),
125 T_(0), tol_(0), m1_(0), m2_(0), m3_(0), nu_(0),
126 ls_maxit_(0), first_print_(true), isConvex_(false),
128 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);
130 state->searchSize = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Initial Trust-Region Parameter", oe3);
131 T_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Maximum Trust-Region Parameter", oe8);
132 tol_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Epsilon Solution Tolerance", oem6);
133 m1_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Upper Threshold for Serious Step", p1);
134 m2_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Lower Threshold for Serious Step", p2);
135 m3_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Upper Threshold for Null Step", p9);
136 nu_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Tolerance for Trust-Region Parameter", oem3);
139 Real coeff = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Distance Measure Coefficient", zero);
140 Real omega = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Locality Measure Coefficient", two);
141 unsigned maxSize = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Maximum Bundle Size", 200);
142 unsigned remSize = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Removal Size for Bundle Update", 2);
143 if ( parlist.sublist(
"Step").sublist(
"Bundle").get(
"Cutting Plane Solver",0) == 1 ) {
145 bundle_ = Teuchos::rcp(
new Bundle_AS<Real>(maxSize,coeff,omega,remSize));
148 bundle_ = Teuchos::rcp(
new Bundle_AS<Real>(maxSize,coeff,omega,remSize));
150 isConvex_ = ((coeff == zero) ?
true :
false);
153 QPtol_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Cutting Plane Tolerance", oem8);
154 QPmaxit_ = parlist.sublist(
"Step").sublist(
"Bundle").get(
"Cutting Plane Iteration Limit", 1000);
158 = parlist.sublist(
"Step").sublist(
"Line Search").get(
"Maximum Number of Function Evaluations",20);
160 lineSearch_ = LineSearchFactory<Real>(parlist);
164 verbosity_ = parlist.sublist(
"General").get(
"Print Verbosity", 0);
172 Real searchSize = state->searchSize;
174 state->searchSize = searchSize;
176 bundle_->initialize(*(state->gradientVec));
180 aggSubGradNew_ = g.
clone();
181 aggSubGradOldNorm_ = algo_state.
gnorm;
184 lineSearch_->initialize(x,x,g,obj,con);
191 first_print_ =
false;
192 QPiter_ = (step_flag_==1 ? 0 :
QPiter_);
193 Real v(0), l(0), u =
T_, gd(0);
194 Real zero(0), two(2), half(0.5);
200 QPiter_ += bundle_->solveDual(state->searchSize,QPmaxit_,QPtol_);
201 bundle_->aggregate(*aggSubGradNew_,aggLinErrNew_,aggDistMeasNew_);
203 if (verbosity_ > 0) {
204 std::cout << std::endl;
205 std::cout <<
" Computation of aggregrate quantities" << std::endl;
207 std::cout <<
" Aggregate linearization error: " << aggLinErrNew_ << std::endl;
208 std::cout <<
" Aggregate distance measure: " << aggDistMeasNew_ << std::endl;
214 s.
set(aggSubGradNew_->dual()); s.
scale(-state->searchSize);
216 if (verbosity_ > 0) {
217 std::cout << std::endl;
218 std::cout <<
" Solve cutting plan subproblem" << std::endl;
219 std::cout <<
" Cutting plan objective value: " << v << std::endl;
220 std::cout <<
" Norm of computed step: " << algo_state.
snorm << std::endl;
221 std::cout <<
" 'Trust-region' radius: " << state->searchSize << std::endl;
231 algo_state.
flag =
true;
235 || std::isnan(aggLinErrNew_)
236 || (std::isnan(aggDistMeasNew_) && !isConvex_)) {
240 algo_state.
flag =
true;
244 y_->set(x); y_->plus(s);
246 valueNew_ = obj.
value(*y_,ftol_);
248 obj.
gradient(*(state->gradientVec),*y_,ftol_);
251 gd = s.
dot(state->gradientVec->dual());
252 linErrNew_ = algo_state.
value - (valueNew_ - gd);
254 Real eps =
static_cast<Real
>(10)*ROL_EPSILON<Real>();
255 Real del = eps*std::max(static_cast<Real>(1),std::abs(algo_state.
value));
256 Real Df = (valueNew_ - algo_state.
value) - del;
259 if (std::abs(Df) < eps && std::abs(Dm) < eps) {
267 bool NS2a = (bundle_->computeAlpha(algo_state.
snorm,linErrNew_) <= m3_*
aggLinErrOld_);
268 bool NS2b = (std::abs(algo_state.
value-valueNew_) <= aggSubGradOldNorm_ +
aggLinErrOld_);
269 if (verbosity_ > 0) {
270 std::cout << std::endl;
271 std::cout <<
" Check for serious/null step" << std::endl;
272 std::cout <<
" Serious step test SS(i): " << SS1 << std::endl;
273 std::cout <<
" -> Left hand side: " << valueNew_-algo_state.
value << std::endl;
274 std::cout <<
" -> Right hand side: " << m1_*v << std::endl;
275 std::cout <<
" Null step test NS(iia): " << NS2a << std::endl;
276 std::cout <<
" -> Left hand side: " << bundle_->computeAlpha(algo_state.
snorm,linErrNew_) << std::endl;
277 std::cout <<
" -> Right hand side: " << m3_*aggLinErrOld_ << std::endl;
278 std::cout <<
" Null step test NS(iib): " << NS2b << std::endl;
279 std::cout <<
" -> Left hand side: " << std::abs(algo_state.
value-valueNew_) << std::endl;
280 std::cout <<
" -> Right hand side: " << aggSubGradOldNorm_ + aggLinErrOld_ << std::endl;
285 bool SS2 = (gd >= m2_*v || state->searchSize >= T_-
nu_);
286 if (verbosity_ > 0) {
287 std::cout <<
" Serious step test SS(iia): " << (gd >= m2_*v) << std::endl;
288 std::cout <<
" -> Left hand side: " << gd << std::endl;
289 std::cout <<
" -> Right hand side: " << m2_*v << std::endl;
290 std::cout <<
" Serious step test SS(iia): " << (state->searchSize >= T_-
nu_) << std::endl;
291 std::cout <<
" -> Left hand side: " << state->searchSize << std::endl;
292 std::cout <<
" -> Right hand side: " << T_-nu_ << std::endl;
297 if (verbosity_ > 0) {
298 std::cout <<
" Serious step taken" << std::endl;
303 l = state->searchSize;
304 state->searchSize = half*(u+l);
305 if (verbosity_ > 0) {
306 std::cout <<
" Increase 'trust-region' radius: " << state->searchSize << std::endl;
311 if ( NS2a || NS2b ) {
315 if (verbosity_ > 0) {
316 std::cout <<
" Null step taken" << std::endl;
321 u = state->searchSize;
322 state->searchSize = half*(u+l);
323 if (verbosity_ > 0) {
324 std::cout <<
" Decrease 'trust-region' radius: " << state->searchSize << std::endl;
331 bool NS3 = (gd - bundle_->computeAlpha(algo_state.
snorm,linErrNew_) >= m2_*v);
332 if (verbosity_ > 0) {
333 std::cout <<
" Null step test NS(iii): " << NS3 << std::endl;
334 std::cout <<
" -> Left hand side: " << gd - bundle_->computeAlpha(algo_state.
snorm,linErrNew_) << std::endl;
335 std::cout <<
" -> Right hand side: " << m2_*v << std::endl;
343 if ( NS2a || NS2b ) {
353 int ls_nfval = 0, ls_ngrad = 0;
354 lineSearch_->run(alpha,valueNew_,ls_nfval,ls_ngrad,gd,s,x,obj,con);
355 if ( ls_nfval == ls_maxit_ ) {
369 u = state->searchSize;
370 state->searchSize = half*(u+l);
375 u = state->searchSize;
376 state->searchSize = half*(u+l);
393 if ( !algo_state.
flag ) {
397 bundle_->reset(*aggSubGradNew_,aggLinErrNew_,algo_state.
snorm);
401 if ( step_flag_==1 ) {
404 Real valueOld = algo_state.
value;
406 bundle_->update(step_flag_,valueNew_-valueOld,algo_state.
snorm,*(state->gradientVec),s);
408 else if ( step_flag_==0 ) {
410 bundle_->update(step_flag_,linErrNew_,algo_state.
snorm,*(state->gradientVec),s);
417 algo_state.
gnorm = (state->gradientVec)->norm();
418 if ( step_flag_==1 ) {
424 std::stringstream hist;
426 hist << std::setw(6) << std::left <<
"iter";
427 hist << std::setw(15) << std::left <<
"value";
428 hist << std::setw(15) << std::left <<
"gnorm";
429 hist << std::setw(15) << std::left <<
"snorm";
430 hist << std::setw(10) << std::left <<
"#fval";
431 hist << std::setw(10) << std::left <<
"#grad";
432 hist << std::setw(15) << std::left <<
"znorm";
433 hist << std::setw(15) << std::left <<
"alpha";
434 hist << std::setw(15) << std::left <<
"TRparam";
435 hist << std::setw(10) << std::left <<
"QPiter";
441 std::stringstream hist;
442 hist <<
"\n" <<
"Bundle Trust-Region Algorithm \n";
448 std::stringstream hist;
449 hist << std::scientific << std::setprecision(6);
450 if ( algo_state.
iter == 0 && first_print_ ) {
452 if ( print_header ) {
456 hist << std::setw(6) << std::left << algo_state.
iter;
457 hist << std::setw(15) << std::left << algo_state.
value;
458 hist << std::setw(15) << std::left << algo_state.
gnorm;
461 if ( step_flag_==1 && algo_state.
iter > 0 ) {
462 if ( print_header ) {
467 hist << std::setw(6) << std::left << algo_state.
iter;
468 hist << std::setw(15) << std::left << algo_state.
value;
469 hist << std::setw(15) << std::left << algo_state.
gnorm;
470 hist << std::setw(15) << std::left << algo_state.
snorm;
471 hist << std::setw(10) << std::left << algo_state.
nfval;
472 hist << std::setw(10) << std::left << algo_state.
ngrad;
475 hist << std::setw(15) << std::left << state->searchSize;
476 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_
Provides the interface for and implements an active set bundle.
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)