44 #ifndef ROL_TRUSTREGIONSTEP_H
45 #define ROL_TRUSTREGIONSTEP_H
128 template <
class Real>
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));
191 ROL::ParameterList &glist = parlist.sublist(
"General");
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));
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));
206 scaleEps_ = glist.get(
"Scale for Epsilon Active Sets", static_cast<Real>(1));
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));
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);
250 Real gtol1 =
scale0_*state->searchSize;
252 Real gtol0 = gtol1 + one;
253 while ( gtol0 > gtol1 ) {
254 obj.
gradient(*(state->gradientVec),x,gtol1);
257 gtol1 =
scale0_*std::min(algo_state.
gnorm,state->searchSize);
262 Real gtol = std::sqrt(ROL_EPSILON<Real>());
263 obj.
gradient(*(state->gradientVec),x,gtol);
290 return xnew_->norm();
330 ROL::ParameterList &glist = parlist.sublist(
"General");
334 secant_ = SecantFactory<Real>(parlist);
363 ROL::ParameterList &glist = parlist.sublist(
"General");
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);
388 Real p1(0.1), oe10(1.e10),
zero(0), one(1), half(0.5), three(3), two(2), six(6);
394 Real htol = std::sqrt(ROL_EPSILON<Real>());
395 Real ftol = p1*ROL_OVERFLOW<Real>();
397 step_state->descentVec = s.
clone();
398 step_state->gradientVec = g.
clone();
415 algo_state.
snorm = oe10;
418 algo_state.
gnorm = ROL_INF<Real>();
425 Ptr<Vector<Real>> v = g.
clone();
426 Ptr<Vector<Real>> hv = x.
clone();
429 catch (std::exception &e) {
435 bool autoRad =
false;
436 if ( step_state->searchSize <=
zero ) {
438 Ptr<Vector<Real>> Bg = g.
clone();
440 secant_->applyB(*Bg,(step_state->gradientVec)->dual());
443 obj.
hessVec(*Bg,(step_state->gradientVec)->dual(),x,htol);
445 Real gBg = Bg->dot(*(step_state->gradientVec));
447 if ( gBg > ROL_EPSILON<Real>() ) {
448 alpha = algo_state.
gnorm*algo_state.
gnorm/gBg;
451 Ptr<Vector<Real>> cp = s.
clone();
452 cp->set((step_state->gradientVec)->dual());
454 Ptr<Vector<Real>> xcp = x.
clone();
461 Real fnew = obj.
value(*xcp,ftol);
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>() ) {
468 step_state->searchSize = std::min(alpha*algo_state.
gnorm,
delMax_);
471 Real b = half*alpha*alpha*gBg;
473 if ( b*b-three*a*c > ROL_EPSILON<Real>() ) {
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 ) {
479 step_state->searchSize = std::min(t1*alpha*algo_state.
gnorm,
delMax_);
483 step_state->searchSize = std::min(t2*alpha*algo_state.
gnorm,
delMax_);
487 step_state->searchSize = std::min(alpha*algo_state.
gnorm,
delMax_);
490 if (step_state->searchSize <= ROL_EPSILON<Real>()*algo_state.
gnorm && autoRad) {
491 step_state->searchSize = one;
498 model_ = makePtr<KelleySachsModel<Real>>(obj,
501 *(step_state->gradientVec),
507 model_ = makePtr<ColemanLiModel<Real>>(obj,
510 *(step_state->gradientVec),
519 model_ = makePtr<LinMoreModel<Real>>(obj,
522 *(step_state->gradientVec),
528 ROL_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
529 ">>> ERROR (TrustRegionStep): Invalid trust-region model!");
533 model_ = makePtr<TrustRegionModel<Real>>(obj,
536 *(step_state->gradientVec),
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);
568 dynamicPtrCast<ColemanLiModel<Real>>(
model_)->setRadius(step_state->searchSize);
603 Real fold = algo_state.
value;
607 s,algo_state.
snorm,fold,*(state->gradientVec),algo_state.
iter,
609 algo_state.
nfval += state->nfval;
610 algo_state.
ngrad += state->ngrad;
611 state->flag =
static_cast<int>(
TRflag_);
620 gp_->set(*(state->gradientVec));
645 algo_state.
value = fnew;
653 std::stringstream hist;
656 hist << std::string(114,
'-') <<
"\n";
658 hist <<
"Trust-Region status output definitions\n\n";
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";
671 hist <<
" tr_flag - Trust-Region flag" <<
"\n";
680 hist <<
" iterCG - Number of Truncated CG iterations\n\n";
681 hist <<
" flagGC - Trust-Region Truncated CG flag" <<
"\n";
688 hist << std::string(114,
'-') <<
"\n";
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";
713 std::stringstream hist;
745 std::stringstream hist;
746 hist << std::scientific << std::setprecision(6);
747 if ( algo_state.
iter == 0 ) {
750 if ( print_header ) {
753 if ( algo_state.
iter == 0 ) {
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;
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_;
std::string ECGFlagToString(ECGFlag cgf)
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...
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.
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)
Defines the linear algebra or vector space interface.
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.
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)
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.
ETrustRegionModel StringToETrustRegionModel(std::string s)
int SPflag_
Subproblem solver termination flag.
ROL::Ptr< StepState< Real > > getState(void)
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.
Real scale0_
Scale for inexact gradient computation.
Ptr< TrustRegionModel< Real > > model_
Container for trust-region model.
ROL::Ptr< Vector< Real > > iterateVec
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.
virtual ~TrustRegionStep()
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 > 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)
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.
Ptr< Vector< Real > > gp_
Container for previous gradient vector.