44#ifndef ROL_TRUSTREGIONALGORITHM_U_DEF_H
45#define ROL_TRUSTREGIONALGORITHM_U_DEF_H
52template<
typename Real>
54 const Ptr<Secant<Real>> &secant )
58 status_->add(makePtr<StatusTest<Real>>(parlist));
61 ParameterList &slist = parlist.sublist(
"Step");
62 ParameterList &trlist = slist.sublist(
"Trust Region");
63 state_->searchSize = trlist.get(
"Initial Radius",
static_cast<Real
>(-1));
64 delMax_ = trlist.get(
"Maximum Radius", ROL_INF<Real>());
65 eta0_ = trlist.get(
"Step Acceptance Threshold",
static_cast<Real
>(0.05));
66 eta1_ = trlist.get(
"Radius Shrinking Threshold",
static_cast<Real
>(0.05));
67 eta2_ = trlist.get(
"Radius Growing Threshold",
static_cast<Real
>(0.9));
68 gamma0_ = trlist.get(
"Radius Shrinking Rate (Negative rho)",
static_cast<Real
>(0.0625));
69 gamma1_ = trlist.get(
"Radius Shrinking Rate (Positive rho)",
static_cast<Real
>(0.25));
70 gamma2_ = trlist.get(
"Radius Growing Rate",
static_cast<Real
>(2.5));
71 TRsafe_ = trlist.get(
"Safeguard Size",
static_cast<Real
>(100.0));
72 eps_ = TRsafe_*ROL_EPSILON<Real>();
74 NMstorage_ = trlist.get(
"Nonmonotone Storage Limit", 0);
75 useNM_ = (NMstorage_ <= 0 ? false :
true);
77 ParameterList &glist = parlist.sublist(
"General");
79 useInexact_.push_back(glist.get(
"Inexact Objective Function",
false));
80 useInexact_.push_back(glist.get(
"Inexact Gradient",
false));
81 useInexact_.push_back(glist.get(
"Inexact Hessian-Times-A-Vector",
false));
83 ParameterList &ilist = trlist.sublist(
"Inexact").sublist(
"Gradient");
84 scale0_ = ilist.get(
"Tolerance Scaling",
static_cast<Real
>(0.1));
85 scale1_ = ilist.get(
"Relative Tolerance",
static_cast<Real
>(2));
87 ParameterList &vlist = trlist.sublist(
"Inexact").sublist(
"Value");
88 scale_ = vlist.get(
"Tolerance Scaling",
static_cast<Real
>(1.e-1));
89 omega_ = vlist.get(
"Exponent",
static_cast<Real
>(0.9));
90 force_ = vlist.get(
"Forcing Sequence Initial Value",
static_cast<Real
>(1.0));
91 updateIter_ = vlist.get(
"Forcing Sequence Update Frequency",
static_cast<int>(10));
92 forceFactor_ = vlist.get(
"Forcing Sequence Reduction Factor",
static_cast<Real
>(0.1));
95 solver_ = TrustRegionUFactory<Real>(parlist);
96 verbosity_ = glist.get(
"Output Level", 0);
98 useSecantPrecond_ = glist.sublist(
"Secant").get(
"Use as Preconditioner",
false);
99 useSecantHessVec_ = glist.sublist(
"Secant").get(
"Use as Hessian",
false);
100 if (secant == nullPtr) {
101 esec_ =
StringToESecant(glist.sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS"));
104 model_ = makePtr<TrustRegionModel_U<Real>>(parlist,secant);
105 printHeader_ = verbosity_ > 2;
108template<
typename Real>
109void TrustRegionAlgorithm<Real>::initialize(
const Vector<Real> &x,
110 const Vector<Real> &g,
112 Objective<Real> &obj,
113 std::ostream &outStream) {
115 Algorithm<Real>::initialize(x,g);
116 solver_->initialize(x,g);
117 model_->initialize(x,g);
119 Real ftol =
static_cast<Real
>(0.1)*ROL_OVERFLOW<Real>();
120 obj.update(x,UpdateType::Initial,state_->iter);
121 state_->value = obj.value(x,ftol);
123 state_->snorm = ROL_INF<Real>();
124 state_->gnorm = ROL_INF<Real>();
125 computeGradient(x,obj);
127 model_->validate(obj,x,g,etr_);
129 if ( state_->searchSize <=
static_cast<Real
>(0) ) {
132 = TRUtils::initialRadius<Real>(nfval,x,*state_->gradientVec,Bg,
133 state_->value,state_->gnorm,obj,*model_,delMax_,
134 outStream,(verbosity_>1));
135 state_->nfval += nfval;
139template<
typename Real>
140Real TrustRegionAlgorithm<Real>::computeValue(
const Vector<Real> &x,
141 Objective<Real> &obj,
144 Real tol(std::sqrt(ROL_EPSILON<Real>())), fval(0);
145 if ( useInexact_[0] ) {
146 if ( !(state_->iter%updateIter_) && (state_->iter != 0) ) {
147 force_ *= forceFactor_;
149 Real eta =
static_cast<Real
>(0.999)*std::min(eta1_,one-eta2_);
150 tol = scale_*std::pow(eta*std::min(pRed,force_),one/omega_);
151 state_->value = obj.value(*state_->iterateVec,tol);
155 obj.update(x,UpdateType::Trial);
156 fval = obj.value(x,tol);
161template<
typename Real>
162void TrustRegionAlgorithm<Real>::computeGradient(
const Vector<Real> &x,
163 Objective<Real> &obj) {
164 if ( useInexact_[1] ) {
166 Real gtol1 = scale0_*state_->searchSize;
167 Real gtol0 = gtol1 + one;
168 while ( gtol0 > gtol1 ) {
169 obj.gradient(*state_->gradientVec,x,gtol1);
170 state_->gnorm = state_->gradientVec->norm();
172 gtol1 = scale0_*std::min(state_->gnorm,state_->searchSize);
176 Real gtol = std::sqrt(ROL_EPSILON<Real>());
177 obj.gradient(*state_->gradientVec,x,gtol);
178 state_->gnorm = state_->gradientVec->norm();
183template<
typename Real>
184void TrustRegionAlgorithm<Real>::run( Vector<Real> &x,
185 const Vector<Real> &g,
186 Objective<Real> &obj,
187 std::ostream &outStream ) {
190 Real ftrial(0), pRed(0), rho(0);
191 Ptr<Vector<Real>> gvec = g.clone();
194 Real rhoNM(0), sigmac(0), sigmar(0);
195 Real fr(state_->value), fc(state_->value), fmin(state_->value);
196 TRUtils::ETRFlag TRflagNM;
200 if (verbosity_ > 0) writeOutput(outStream,
true);
202 while (status_->check(*state_)) {
204 model_->setData(obj,x,*state_->gradientVec);
207 SPflag_ = 0; SPiter_ = 0;
208 solver_->solve(*state_->stepVec,state_->snorm,pRed,SPflag_,SPiter_,
209 state_->searchSize,*model_);
211 x.plus(*state_->stepVec);
212 ftrial = computeValue(x,obj,pRed);
214 TRflag_ = TRUtils::SUCCESS;
215 TRUtils::analyzeRatio<Real>(rho,TRflag_,state_->value,ftrial,pRed,eps_,outStream,verbosity_>1);
217 TRUtils::analyzeRatio<Real>(rhoNM,TRflagNM,fr,ftrial,pRed+sigmar,eps_,outStream,verbosity_>1);
218 TRflag_ = (rho < rhoNM ? TRflagNM : TRflag_);
219 rho = (rho < rhoNM ? rhoNM : rho );
224 if ((rho < eta0_ && TRflag_ == TRUtils::SUCCESS)
226 x.set(*state_->iterateVec);
227 obj.update(x,UpdateType::Revert,state_->iter);
228 if (rho <
zero && TRflag_ != TRUtils::TRNAN) {
230 state_->searchSize = TRUtils::interpolateRadius<Real>(*state_->gradientVec,*state_->stepVec,
231 state_->snorm,pRed,state_->value,ftrial,state_->searchSize,gamma0_,gamma1_,eta2_,
232 outStream,verbosity_>1);
235 state_->searchSize = gamma1_*std::min(state_->snorm,state_->searchSize);
237 if (useInexact_[1]) computeGradient(x,obj);
239 else if ((rho >= eta0_ && TRflag_ != TRUtils::NPOSPREDNEG)
240 || (TRflag_ == TRUtils::POSPREDNEG)) {
241 state_->iterateVec->set(x);
242 state_->value = ftrial;
243 obj.update(x,UpdateType::Accept,state_->iter);
245 sigmac += pRed; sigmar += pRed;
246 if (ftrial < fmin) { fmin = ftrial; fc = fmin; sigmac =
zero; L = 0; }
249 if (ftrial > fc) { fc = ftrial; sigmac =
zero; }
250 if (L == NMstorage_) { fr = fc; sigmar = sigmac; }
254 if (rho >= eta2_) state_->searchSize = std::min(gamma2_*state_->searchSize, delMax_);
256 gvec->set(*state_->gradientVec);
257 computeGradient(x,obj);
259 model_->update(x,*state_->stepVec,*gvec,*state_->gradientVec,
260 state_->snorm,state_->iter);
263 if (verbosity_ > 0) writeOutput(outStream,printHeader_);
265 if (verbosity_ > 0) Algorithm<Real>::writeExitStatus(outStream);
268template<
typename Real>
269void TrustRegionAlgorithm<Real>::writeHeader( std::ostream& os )
const {
270 std::stringstream hist;
272 hist << std::string(114,
'-') << std::endl;
273 hist <<
"Trust-Region status output definitions" << std::endl << std::endl;
274 hist <<
" iter - Number of iterates (steps taken)" << std::endl;
275 hist <<
" value - Objective function value" << std::endl;
276 hist <<
" gnorm - Norm of the gradient" << std::endl;
277 hist <<
" snorm - Norm of the step (update to optimization vector)" << std::endl;
278 hist <<
" delta - Trust-Region radius" << std::endl;
279 hist <<
" #fval - Number of times the objective function was evaluated" << std::endl;
280 hist <<
" #grad - Number of times the gradient was computed" << std::endl;
282 hist <<
" tr_flag - Trust-Region flag" << std::endl;
283 for(
int flag = TRUtils::SUCCESS; flag != TRUtils::UNDEFINED; ++flag ) {
285 << TRUtils::ETRFlagToString(
static_cast<TRUtils::ETRFlag
>(flag)) << std::endl;
287 if( etr_ == TRUSTREGION_U_TRUNCATEDCG ) {
289 hist <<
" iterCG - Number of Truncated CG iterations" << std::endl << std::endl;
290 hist <<
" flagGC - Trust-Region Truncated CG flag" << std::endl;
296 else if( etr_ == TRUSTREGION_U_SPG ) {
298 hist <<
" iterCG - Number of spectral projected gradient iterations" << std::endl << std::endl;
299 hist <<
" flagGC - Trust-Region spectral projected gradient flag" << std::endl;
301 hist << std::string(114,
'-') << std::endl;
304 hist << std::setw(6) << std::left <<
"iter";
305 hist << std::setw(15) << std::left <<
"value";
306 hist << std::setw(15) << std::left <<
"gnorm";
307 hist << std::setw(15) << std::left <<
"snorm";
308 hist << std::setw(15) << std::left <<
"delta";
309 hist << std::setw(10) << std::left <<
"#fval";
310 hist << std::setw(10) << std::left <<
"#grad";
311 hist << std::setw(10) << std::left <<
"tr_flag";
312 if ( etr_ == TRUSTREGION_U_TRUNCATEDCG ) {
313 hist << std::setw(10) << std::left <<
"iterCG";
314 hist << std::setw(10) << std::left <<
"flagCG";
316 else if (etr_ == TRUSTREGION_U_SPG) {
317 hist << std::setw(10) << std::left <<
"iterSPG";
318 hist << std::setw(10) << std::left <<
"flagSPG";
324template<
typename Real>
325void TrustRegionAlgorithm<Real>::writeName( std::ostream& os )
const {
326 std::stringstream hist;
328 if ( useSecantPrecond_ || useSecantHessVec_ ) {
329 if ( useSecantPrecond_ && !useSecantHessVec_ ) {
330 hist <<
" with " <<
ESecantToString(esec_) <<
" Preconditioning" << std::endl;
332 else if ( !useSecantPrecond_ && useSecantHessVec_ ) {
333 hist <<
" with " <<
ESecantToString(esec_) <<
" Hessian Approximation" << std::endl;
336 hist <<
" with " <<
ESecantToString(esec_) <<
" Preconditioning and Hessian Approximation" << std::endl;
345template<
typename Real>
346void TrustRegionAlgorithm<Real>::writeOutput(std::ostream& os,
bool print_header)
const {
347 std::stringstream hist;
348 hist << std::scientific << std::setprecision(6);
349 if ( state_->iter == 0 ) {
352 if ( print_header ) {
355 if ( state_->iter == 0 ) {
357 hist << std::setw(6) << std::left << state_->iter;
358 hist << std::setw(15) << std::left << state_->value;
359 hist << std::setw(15) << std::left << state_->gnorm;
360 hist << std::setw(15) << std::left <<
"---";
361 hist << std::setw(15) << std::left << state_->searchSize;
362 hist << std::setw(10) << std::left << state_->nfval;
363 hist << std::setw(10) << std::left << state_->ngrad;
364 hist << std::setw(10) << std::left <<
"---";
365 if ( etr_ == TRUSTREGION_U_TRUNCATEDCG || etr_ == TRUSTREGION_U_SPG ) {
366 hist << std::setw(10) << std::left <<
"---";
367 hist << std::setw(10) << std::left <<
"---";
373 hist << std::setw(6) << std::left << state_->iter;
374 hist << std::setw(15) << std::left << state_->value;
375 hist << std::setw(15) << std::left << state_->gnorm;
376 hist << std::setw(15) << std::left << state_->snorm;
377 hist << std::setw(15) << std::left << state_->searchSize;
378 hist << std::setw(10) << std::left << state_->nfval;
379 hist << std::setw(10) << std::left << state_->ngrad;
380 hist << std::setw(10) << std::left << TRflag_;
381 if ( etr_ == TRUSTREGION_U_TRUNCATEDCG || etr_ == TRUSTREGION_U_SPG ) {
382 hist << std::setw(10) << std::left << SPiter_;
383 hist << std::setw(10) << std::left << SPflag_;
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0 zero)()
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.
Contains definitions of enums for trust region algorithms.
Provides an interface to run unconstrained optimization algorithms.
TrustRegionAlgorithm(ParameterList &parlist, const Ptr< Secant< Real > > &secant=nullPtr)
std::string NumberToString(T Number)
ETrustRegionU StringToETrustRegionU(std::string s)
ESecant StringToESecant(std::string s)
std::string ETrustRegionUToString(ETrustRegionU tr)
ECGFlag
Enumation of flags used by conjugate gradient methods.
std::string ESecantToString(ESecant tr)
std::string ECGFlagToString(ECGFlag cgf)