44 #ifndef ROL_TRUSTREGION_H 45 #define ROL_TRUSTREGION_H 89 : ftol_old_(
ROL_OVERFLOW<Real>()), cnt_(0), verbosity_(0) {
91 Teuchos::ParameterList list = parlist.sublist(
"Step").sublist(
"Trust Region");
93 delmax_ = list.get(
"Maximum Radius", static_cast<Real>(5000.0));
94 eta0_ = list.get(
"Step Acceptance Threshold", static_cast<Real>(0.05));
95 eta1_ = list.get(
"Radius Shrinking Threshold", static_cast<Real>(0.05));
96 eta2_ = list.get(
"Radius Growing Threshold", static_cast<Real>(0.9));
97 gamma0_ = list.get(
"Radius Shrinking Rate (Negative rho)", static_cast<Real>(0.0625));
98 gamma1_ = list.get(
"Radius Shrinking Rate (Positive rho)", static_cast<Real>(0.25));
99 gamma2_ = list.get(
"Radius Growing Rate", static_cast<Real>(2.5));
100 mu0_ = list.get(
"Sufficient Decrease Parameter", static_cast<Real>(1.e-4));
101 TRsafe_ = list.get(
"Safeguard Size", static_cast<Real>(100.0));
102 eps_ = TRsafe_*ROL_EPSILON<Real>();
104 Teuchos::ParameterList &glist = parlist.sublist(
"General");
106 useInexact_.push_back(glist.get(
"Inexact Objective Function",
false));
107 useInexact_.push_back(glist.get(
"Inexact Gradient",
false));
108 useInexact_.push_back(glist.get(
"Inexact Hessian-Times-A-Vector",
false));
110 Teuchos::ParameterList &ilist = list.sublist(
"Inexact").sublist(
"Value");
111 scale_ = ilist.get(
"Tolerance Scaling", static_cast<Real>(1.e-1));
112 omega_ = ilist.get(
"Exponent", static_cast<Real>(0.9));
113 force_ = ilist.get(
"Forcing Sequence Initial Value", static_cast<Real>(1.0));
114 updateIter_ = ilist.get(
"Forcing Sequence Update Frequency", static_cast<int>(10));
115 forceFactor_ = ilist.get(
"Forcing Sequence Reduction Factor", static_cast<Real>(0.1));
117 verbosity_ = glist.get(
"Print Verbosity", 0);
139 Real tol = std::sqrt(ROL_EPSILON<Real>());
140 const Real one(1), oe4(1.e4), zero(0);
146 Real fold1 = fold, ftol = tol, TOL(1.e-2);
147 if ( useInexact_[0] ) {
148 if ( !(cnt_%updateIter_) && (cnt_ != 0) ) {
151 Real c = scale_*std::max(TOL,std::min(one,oe4*std::max(pRed_,std::sqrt(ROL_EPSILON<Real>()))));
152 ftol = c*std::pow(std::min(eta1_,one-eta2_)
153 *std::min(std::max(pRed_,std::sqrt(ROL_EPSILON<Real>())),force_),one/omega_);
154 if ( ftol_old_ > ftol || cnt_ == 0 ) {
156 fold1 = obj.
value(x,ftol_old_);
161 prim_->set(x); prim_->plus(s);
163 fnew = obj.
value(*prim_,ftol);
166 Real aRed = fold1 - fnew;
178 if ( verbosity_ > 0 ) {
179 std::cout << std::endl;
180 std::cout <<
" Computation of actual and predicted reduction" << std::endl;
181 std::cout <<
" Current objective function value: " << fold1 << std::endl;
182 std::cout <<
" New objective function value: " << fnew << std::endl;
183 std::cout <<
" Actual reduction: " << aRed << std::endl;
184 std::cout <<
" Predicted reduction: " << pRed_ << std::endl;
188 Real EPS = eps_*((one > std::abs(fold1)) ? one : std::abs(fold1));
189 Real aRed_safe = aRed + EPS, pRed_safe = pRed_ + EPS;
191 if (((std::abs(aRed_safe) < eps_) && (std::abs(pRed_safe) < eps_)) || aRed == pRed_) {
195 else if ( std::isnan(aRed_safe) || std::isnan(pRed_safe) ) {
200 rho = aRed_safe/pRed_safe;
201 if (pRed_safe < zero && aRed_safe > zero) {
204 else if (aRed_safe <= zero && pRed_safe > zero) {
207 else if (aRed_safe <= zero && pRed_safe < zero) {
215 if ( verbosity_ > 0 ) {
216 std::cout <<
" Safeguard: " << eps_ << std::endl;
217 std::cout <<
" Actual reduction with safeguard: " << aRed_safe << std::endl;
218 std::cout <<
" Predicted reduction with safeguard: " << pRed_safe << std::endl;
219 std::cout <<
" Ratio of actual and predicted reduction: " << rho << std::endl;
220 std::cout <<
" Trust-region flag: " << flagTR << std::endl;
232 if ( rho >= eta0_ && (std::abs(aRed_safe) > eps_) ) {
235 prim_->axpy(-one,g.
dual());
239 Real pgnorm = prim_->norm();
241 prim_->set(g.
dual());
243 Real lam = std::min(one, del/prim_->norm());
249 pgnorm *= prim_->norm();
251 decr = ( aRed_safe >= mu0_*eta0_*pgnorm );
254 if ( verbosity_ > 0 ) {
255 std::cout <<
" Decrease lower bound (constraints): " << 0.1*eta0_*pgnorm << std::endl;
256 std::cout <<
" Trust-region flag (constraints): " << flagTR << std::endl;
257 std::cout <<
" Is step feasible: " << bnd.
isFeasible(x) << std::endl;
268 if ( verbosity_ > 0 ) {
269 std::cout <<
" Norm of step: " << snorm << std::endl;
270 std::cout <<
" Trust-region radius before update: " << del << std::endl;
278 gs = dual_->dot(s.
dual());
283 Real modelVal = model.
value(s,tol);
285 Real theta = (one-
eta2_)*gs/((one-eta2_)*(fold1+gs)+eta2_*modelVal-fnew);
286 del = std::min(gamma1_*std::min(snorm,del),std::max(gamma0_,theta)*del);
289 del = gamma1_*std::min(snorm,del);
298 del = std::min(gamma2_*del,delmax_);
301 if ( verbosity_ > 0 ) {
302 std::cout <<
" Trust-region radius after update: " << del << std::endl;
303 std::cout << std::endl;
Provides the interface to evaluate objective functions.
TrustRegion(Teuchos::ParameterList &parlist)
virtual void plus(const Vector &x)=0
Compute , where .
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
virtual void initialize(const Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g)
Contains definitions of custom data types in ROL.
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
Provides interface for and implements trust-region subproblem solvers.
Provides the interface to evaluate trust-region model functions.
Defines the linear algebra or vector space interface.
virtual void updatePredictedReduction(Real &pred, const Vector< Real > &s)
virtual Real dot(const Vector &x) const =0
Compute where .
virtual void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=0)
Set variables to zero if they correspond to the -active set.
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
Teuchos::RCP< Vector< Real > > prim_
bool isActivated(void)
Check if bounds are on.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
ETrustRegionModel StringToETrustRegionModel(std::string s)
Real getPredictedReduction(void) const
virtual Real value(const Vector< Real > &s, Real &tol)
Compute value.
virtual const Teuchos::RCP< const Vector< Real > > getGradient(void) const
void setPredictedReduction(const Real pRed)
Real ROL_OVERFLOW(void)
Platform-dependent maximum double.
ETrustRegionModel
Enumeration of trust-region model types.
Provides the interface to apply upper and lower bound constraints.
Teuchos::RCP< Vector< Real > > dual_
virtual void updateActualReduction(Real &ared, const Vector< Real > &s)
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
ETrustRegionFlag
Enumation of flags used by trust-region solvers.
ETrustRegionModel TRmodel_
virtual bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
Contains definitions of enums for trust region algorithms.
std::vector< bool > useInexact_
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
virtual void run(Vector< Real > &s, Real &snorm, int &iflag, int &iter, const Real del, TrustRegionModel< Real > &model)=0
virtual void update(Vector< Real > &x, Real &fnew, Real &del, int &nfval, int &ngrad, ETrustRegionFlag &flagTR, const Vector< Real > &s, const Real snorm, const Real fold, const Vector< Real > &g, int iter, Objective< Real > &obj, BoundConstraint< Real > &bnd, TrustRegionModel< Real > &model)