44 #ifndef ROL_TRUSTREGION_H
45 #define ROL_TRUSTREGION_H
96 ROL::ParameterList list = parlist.sublist(
"Step").sublist(
"Trust Region");
98 eta0_ = list.get(
"Step Acceptance Threshold", static_cast<Real>(0.05));
99 eta1_ = list.get(
"Radius Shrinking Threshold", static_cast<Real>(0.05));
100 eta2_ = list.get(
"Radius Growing Threshold", static_cast<Real>(0.9));
101 gamma0_ = list.get(
"Radius Shrinking Rate (Negative rho)", static_cast<Real>(0.0625));
102 gamma1_ = list.get(
"Radius Shrinking Rate (Positive rho)", static_cast<Real>(0.25));
103 gamma2_ = list.get(
"Radius Growing Rate", static_cast<Real>(2.5));
104 mu0_ = list.get(
"Sufficient Decrease Parameter", static_cast<Real>(1.e-4));
105 TRsafe_ = list.get(
"Safeguard Size", static_cast<Real>(100.0));
108 ROL::ParameterList &glist = parlist.sublist(
"General");
110 useInexact_.push_back(glist.get(
"Inexact Objective Function",
false));
111 useInexact_.push_back(glist.get(
"Inexact Gradient",
false));
112 useInexact_.push_back(glist.get(
"Inexact Hessian-Times-A-Vector",
false));
114 ROL::ParameterList &ilist = list.sublist(
"Inexact").sublist(
"Value");
115 scale_ = ilist.get(
"Tolerance Scaling", static_cast<Real>(1.e-1));
116 omega_ = ilist.get(
"Exponent", static_cast<Real>(0.9));
117 force_ = ilist.get(
"Forcing Sequence Initial Value", static_cast<Real>(1.0));
118 updateIter_ = ilist.get(
"Forcing Sequence Update Frequency", static_cast<int>(10));
119 forceFactor_ = ilist.get(
"Forcing Sequence Reduction Factor", static_cast<Real>(0.1));
123 max_fval_ = list.sublist(
"Post-Smoothing").get(
"Function Evaluation Limit", 20);
124 alpha_init_ = list.sublist(
"Post-Smoothing").get(
"Initial Step Size", static_cast<Real>(1));
125 mu_ = list.sublist(
"Post-Smoothing").get(
"Tolerance", static_cast<Real>(0.9999));
126 beta_ = list.sublist(
"Post-Smoothing").get(
"Rate", static_cast<Real>(0.01));
149 Real tol = std::sqrt(ROL_EPSILON<Real>());
150 const Real one(1),
zero(0);
156 Real fold1 = fold, ftol = tol;
170 Real eta =
static_cast<Real
>(0.999)*std::min(
eta1_,one-
eta2_);
185 Real aRed = fold1 - fnew;
198 std::cout << std::endl;
199 std::cout <<
" Computation of actual and predicted reduction" << std::endl;
200 std::cout <<
" Current objective function value: " << fold1 << std::endl;
201 std::cout <<
" New objective function value: " << fnew << std::endl;
202 std::cout <<
" Actual reduction: " << aRed << std::endl;
203 std::cout <<
" Predicted reduction: " <<
pRed_ << std::endl;
207 Real EPS =
eps_*((one > std::abs(fold1)) ? one : std::abs(fold1));
208 Real aRed_safe = aRed + EPS, pRed_safe =
pRed_ + EPS;
210 if (((std::abs(aRed_safe) <
eps_) && (std::abs(pRed_safe) <
eps_)) || aRed ==
pRed_) {
214 else if ( std::isnan(aRed_safe) || std::isnan(pRed_safe) ) {
219 rho = aRed_safe/pRed_safe;
220 if (pRed_safe < zero && aRed_safe > zero) {
223 else if (aRed_safe <= zero && pRed_safe > zero) {
226 else if (aRed_safe <= zero && pRed_safe < zero) {
235 std::cout <<
" Safeguard: " <<
eps_ << std::endl;
236 std::cout <<
" Actual reduction with safeguard: " << aRed_safe << std::endl;
237 std::cout <<
" Predicted reduction with safeguard: " << pRed_safe << std::endl;
238 std::cout <<
" Ratio of actual and predicted reduction: " << rho << std::endl;
239 std::cout <<
" Trust-region flag: " << flagTR << std::endl;
250 if ( rho >=
eta0_ && (std::abs(aRed_safe) >
eps_) ) {
257 Real pgnorm =
prim_->norm();
261 Real lam = std::min(one, del/
prim_->norm());
267 pgnorm *=
prim_->norm();
269 decr = ( aRed_safe >=
mu0_*pgnorm );
273 std::cout <<
" Decrease lower bound (constraints): " <<
mu0_*pgnorm << std::endl;
274 std::cout <<
" Trust-region flag (constraints): " << flagTR << std::endl;
275 std::cout <<
" Is step feasible: " << bnd.
isFeasible(x) << std::endl;
287 std::cout <<
" Norm of step: " << snorm << std::endl;
288 std::cout <<
" Trust-region radius before update: " << del << std::endl;
301 Real modelVal = model.
value(s,tol);
303 Real theta = (one-
eta2_)*gs/((one-
eta2_)*(fold1+gs)+
eta2_*modelVal-fnew);
304 del = std::min(
gamma1_*std::min(snorm,del),std::max(
gamma0_,theta)*del);
306 std::cout <<
" Interpolation model value: " << modelVal << std::endl;
307 std::cout <<
" Interpolation step length: " << theta << std::endl;
311 del =
gamma1_*std::min(snorm,del);
336 while ( (ftmp-fnew) >=
mu_*aRed ) {
350 if (std::isnan(ftmp)) {
352 del =
gamma1_*std::min(snorm,del);
353 rho =
static_cast<Real
>(-1);
373 std::cout <<
" Trust-region radius after update: " << del << std::endl;
374 std::cout << std::endl;
Provides the interface to evaluate objective functions.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void plus(const Vector &x)=0
Compute , where .
bool isActivated(void) const
Check if bounds are on.
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.
Provides interface for and implements trust-region subproblem solvers.
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=0)
Set variables to zero if they correspond to the -active set.
Provides the interface to evaluate trust-region model functions.
Defines the linear algebra or vector space interface.
virtual const Ptr< const Vector< Real > > getGradient(void) const
virtual void updatePredictedReduction(Real &pred, const Vector< Real > &s)
Real alpha_init_
Initial line-search parameter for projected methods.
ROL::Ptr< Vector< Real > > xtmp_
virtual Real dot(const Vector &x) const =0
Compute where .
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
virtual void dualTransform(Vector< Real > &tv, const Vector< Real > &v)
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Real getPredictedReduction(void) const
ETrustRegionModel StringToETrustRegionModel(std::string s)
Real mu_
Post-Smoothing tolerance for projected methods.
Real beta_
Post-Smoothing rate for projected methods.
virtual Real value(const Vector< Real > &s, Real &tol)
Compute value.
int max_fval_
Maximum function evaluations in line-search for projected methods.
void setPredictedReduction(const Real pRed)
Real ROL_OVERFLOW(void)
Platform-dependent maximum double.
ETrustRegionModel
Enumeration of trust-region model types.
TrustRegion(ROL::ParameterList &parlist)
Provides the interface to apply upper and lower bound constraints.
virtual void set(const Vector &x)
Set where .
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_
ROL::Ptr< Vector< Real > > dual_
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)
ROL::Ptr< Vector< Real > > prim_