10#ifndef ROL_TRUSTREGIONALGORITHM_U_DEF_H
11#define ROL_TRUSTREGIONALGORITHM_U_DEF_H
18template<
typename Real>
20 const Ptr<Secant<Real>> &secant )
24 status_->add(makePtr<StatusTest<Real>>(parlist));
27 ParameterList &slist = parlist.sublist(
"Step");
28 ParameterList &trlist = slist.sublist(
"Trust Region");
29 state_->searchSize = trlist.get(
"Initial Radius",
static_cast<Real
>(-1));
30 delMax_ = trlist.get(
"Maximum Radius", ROL_INF<Real>());
31 eta0_ = trlist.get(
"Step Acceptance Threshold",
static_cast<Real
>(0.05));
32 eta1_ = trlist.get(
"Radius Shrinking Threshold",
static_cast<Real
>(0.05));
33 eta2_ = trlist.get(
"Radius Growing Threshold",
static_cast<Real
>(0.9));
34 gamma0_ = trlist.get(
"Radius Shrinking Rate (Negative rho)",
static_cast<Real
>(0.0625));
35 gamma1_ = trlist.get(
"Radius Shrinking Rate (Positive rho)",
static_cast<Real
>(0.25));
36 gamma2_ = trlist.get(
"Radius Growing Rate",
static_cast<Real
>(2.5));
37 TRsafe_ = trlist.get(
"Safeguard Size",
static_cast<Real
>(100.0));
38 eps_ = TRsafe_*ROL_EPSILON<Real>();
40 NMstorage_ = trlist.get(
"Nonmonotone Storage Limit", 0);
41 useNM_ = (NMstorage_ <= 0 ? false :
true);
43 ParameterList &glist = parlist.sublist(
"General");
45 useInexact_.push_back(glist.get(
"Inexact Objective Function",
false));
46 useInexact_.push_back(glist.get(
"Inexact Gradient",
false));
47 useInexact_.push_back(glist.get(
"Inexact Hessian-Times-A-Vector",
false));
49 ParameterList &ilist = trlist.sublist(
"Inexact").sublist(
"Gradient");
50 scale0_ = ilist.get(
"Tolerance Scaling",
static_cast<Real
>(0.1));
51 scale1_ = ilist.get(
"Relative Tolerance",
static_cast<Real
>(2));
53 ParameterList &vlist = trlist.sublist(
"Inexact").sublist(
"Value");
54 scale_ = vlist.get(
"Tolerance Scaling",
static_cast<Real
>(1.e-1));
55 omega_ = vlist.get(
"Exponent",
static_cast<Real
>(0.9));
56 force_ = vlist.get(
"Forcing Sequence Initial Value",
static_cast<Real
>(1.0));
57 updateIter_ = vlist.get(
"Forcing Sequence Update Frequency",
static_cast<int>(10));
58 forceFactor_ = vlist.get(
"Forcing Sequence Reduction Factor",
static_cast<Real
>(0.1));
60 std::string solverName = trlist.get(
"Subproblem Solver",
"Dogleg");
62 solver_ = TrustRegionUFactory<Real>(parlist);
63 verbosity_ = glist.get(
"Output Level", 0);
65 useSecantPrecond_ = glist.sublist(
"Secant").get(
"Use as Preconditioner",
false);
66 useSecantHessVec_ = glist.sublist(
"Secant").get(
"Use as Hessian",
false);
67 if (secant == nullPtr) {
68 std::string secantType = glist.sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS");
72 model_ = makePtr<TrustRegionModel_U<Real>>(parlist,secant);
73 printHeader_ = verbosity_ > 2;
76template<
typename Real>
77void TrustRegionAlgorithm<Real>::initialize(
const Vector<Real> &x,
78 const Vector<Real> &g,
81 std::ostream &outStream) {
84 Algorithm<Real>::initialize(x,g);
85 solver_->initialize(x,g);
86 model_->initialize(x,g);
88 Real ftol =
static_cast<Real
>(0.1)*ROL_OVERFLOW<Real>();
89 obj.update(x,UpdateType::Initial,state_->iter);
90 state_->value = obj.value(x,ftol);
92 state_->snorm = ROL_INF<Real>();
93 state_->gnorm = ROL_INF<Real>();
94 Real Delta = state_->searchSize;
95 if (Delta <=
zero) state_->searchSize = 1e2*x.norm();
96 computeGradient(x,obj,
true);
98 model_->validate(obj,x,g,etr_);
100 if ( Delta <=
zero ) {
103 = TRUtils::initialRadius<Real>(nfval,x,*state_->gradientVec,Bg,
104 state_->value,state_->gnorm,gtol_,obj,*model_,delMax_,
105 outStream,(verbosity_>1));
106 state_->nfval += nfval;
110template<
typename Real>
111Real TrustRegionAlgorithm<Real>::computeValue(
const Vector<Real> &x,
112 Objective<Real> &obj,
115 Real tol(std::sqrt(ROL_EPSILON<Real>())), fval(0);
116 if ( useInexact_[0] ) {
117 if ( !(state_->iter%updateIter_) && (state_->iter != 0) ) {
118 force_ *= forceFactor_;
120 Real eta =
static_cast<Real
>(0.999)*std::min(eta1_,one-eta2_);
121 tol = scale_*std::pow(eta*std::min(pRed,force_),one/omega_);
122 state_->value = obj.value(*state_->iterateVec,tol);
126 obj.update(x,UpdateType::Trial);
127 fval = obj.value(x,tol);
132template<
typename Real>
133void TrustRegionAlgorithm<Real>::computeGradient(
const Vector<Real> &x,
134 Objective<Real> &obj,
136 if ( useInexact_[1] ) {
137 Real gtol0 = scale0_*state_->searchSize;
138 if (accept) gtol_ = gtol0 +
static_cast<Real
>(1);
139 else gtol0 = scale0_*std::min(state_->gnorm,state_->searchSize);
140 while ( gtol_ > gtol0 ) {
142 obj.gradient(*state_->gradientVec,x,gtol_); state_->ngrad++;
143 state_->gnorm = state_->gradientVec->norm();
144 gtol0 = scale0_*std::min(state_->gnorm,state_->searchSize);
149 gtol_ = std::sqrt(ROL_EPSILON<Real>());
150 obj.gradient(*state_->gradientVec,x,gtol_); state_->ngrad++;
151 state_->gnorm = state_->gradientVec->norm();
156template<
typename Real>
157void TrustRegionAlgorithm<Real>::run( Vector<Real> &x,
158 const Vector<Real> &g,
159 Objective<Real> &obj,
160 std::ostream &outStream ) {
163 Real ftrial(0), pRed(0), rho(0);
164 Ptr<Vector<Real>> gvec = g.clone();
167 Real rhoNM(0), sigmac(0), sigmar(0);
168 Real fr(state_->value), fc(state_->value), fmin(state_->value);
169 TRUtils::ETRFlag TRflagNM;
173 if (verbosity_ > 0) writeOutput(outStream,
true);
175 while (status_->check(*state_)) {
177 model_->setData(obj,x,*state_->gradientVec,gtol_);
180 SPflag_ = 0; SPiter_ = 0;
181 solver_->solve(*state_->stepVec,state_->snorm,pRed,SPflag_,SPiter_,
182 state_->searchSize,*model_);
184 x.plus(*state_->stepVec);
185 ftrial = computeValue(x,obj,pRed);
187 TRflag_ = TRUtils::SUCCESS;
188 TRUtils::analyzeRatio<Real>(rho,TRflag_,state_->value,ftrial,pRed,eps_,outStream,verbosity_>1);
190 TRUtils::analyzeRatio<Real>(rhoNM,TRflagNM,fr,ftrial,pRed+sigmar,eps_,outStream,verbosity_>1);
191 TRflag_ = (rho < rhoNM ? TRflagNM : TRflag_);
192 rho = (rho < rhoNM ? rhoNM : rho );
197 if ((rho < eta0_ && TRflag_ == TRUtils::SUCCESS)
199 x.set(*state_->iterateVec);
200 obj.update(x,UpdateType::Revert,state_->iter);
201 if (rho <
zero && TRflag_ != TRUtils::TRNAN) {
203 state_->searchSize = TRUtils::interpolateRadius<Real>(*state_->gradientVec,*state_->stepVec,
204 state_->snorm,pRed,state_->value,ftrial,state_->searchSize,gamma0_,gamma1_,eta2_,
205 outStream,verbosity_>1);
208 state_->searchSize = gamma1_*std::min(state_->snorm,state_->searchSize);
210 computeGradient(x,obj,
false);
212 else if ((rho >= eta0_ && TRflag_ != TRUtils::NPOSPREDNEG)
213 || (TRflag_ == TRUtils::POSPREDNEG)) {
214 state_->iterateVec->set(x);
215 state_->value = ftrial;
216 obj.update(x,UpdateType::Accept,state_->iter);
218 sigmac += pRed; sigmar += pRed;
219 if (ftrial < fmin) { fmin = ftrial; fc = fmin; sigmac =
zero; L = 0; }
222 if (ftrial > fc) { fc = ftrial; sigmac =
zero; }
223 if (L == NMstorage_) { fr = fc; sigmar = sigmac; }
227 if (rho >= eta2_) state_->searchSize = std::min(gamma2_*state_->searchSize, delMax_);
229 gvec->set(*state_->gradientVec);
230 computeGradient(x,obj,
true);
232 model_->update(x,*state_->stepVec,*gvec,*state_->gradientVec,
233 state_->snorm,state_->iter);
236 if (verbosity_ > 0) writeOutput(outStream,printHeader_);
238 if (verbosity_ > 0) Algorithm<Real>::writeExitStatus(outStream);
241template<
typename Real>
242void TrustRegionAlgorithm<Real>::writeHeader( std::ostream& os )
const {
243 std::ios_base::fmtflags osFlags(os.flags());
245 os << std::string(114,
'-') << std::endl;
246 os <<
"Trust-Region status output definitions" << std::endl << std::endl;
247 os <<
" iter - Number of iterates (steps taken)" << std::endl;
248 os <<
" value - Objective function value" << std::endl;
249 os <<
" gnorm - Norm of the gradient" << std::endl;
250 os <<
" snorm - Norm of the step (update to optimization vector)" << std::endl;
251 os <<
" delta - Trust-Region radius" << std::endl;
252 os <<
" #fval - Number of times the objective function was evaluated" << std::endl;
253 os <<
" #grad - Number of times the gradient was computed" << std::endl;
255 os <<
" tr_flag - Trust-Region flag" << std::endl;
256 for(
int flag = TRUtils::SUCCESS; flag != TRUtils::UNDEFINED; ++flag ) {
258 << TRUtils::ETRFlagToString(
static_cast<TRUtils::ETRFlag
>(flag)) << std::endl;
260 if( etr_ == TRUSTREGION_U_TRUNCATEDCG ) {
262 os <<
" iterCG - Number of Truncated CG iterations" << std::endl << std::endl;
263 os <<
" flagGC - Trust-Region Truncated CG flag" << std::endl;
269 else if( etr_ == TRUSTREGION_U_SPG ) {
271 os <<
" iterCG - Number of spectral projected gradient iterations" << std::endl << std::endl;
272 os <<
" flagGC - Trust-Region spectral projected gradient flag" << std::endl;
274 os << std::string(114,
'-') << std::endl;
277 os << std::setw(6) << std::left <<
"iter";
278 os << std::setw(15) << std::left <<
"value";
279 os << std::setw(15) << std::left <<
"gnorm";
280 os << std::setw(15) << std::left <<
"snorm";
281 os << std::setw(15) << std::left <<
"delta";
282 os << std::setw(10) << std::left <<
"#fval";
283 os << std::setw(10) << std::left <<
"#grad";
284 os << std::setw(10) << std::left <<
"tr_flag";
285 if ( etr_ == TRUSTREGION_U_TRUNCATEDCG ) {
286 os << std::setw(10) << std::left <<
"iterCG";
287 os << std::setw(10) << std::left <<
"flagCG";
289 else if (etr_ == TRUSTREGION_U_SPG) {
290 os << std::setw(10) << std::left <<
"iterSPG";
291 os << std::setw(10) << std::left <<
"flagSPG";
297template<
typename Real>
298void TrustRegionAlgorithm<Real>::writeName( std::ostream& os )
const {
299 std::ios_base::fmtflags osFlags(os.flags());
301 if ( useSecantPrecond_ || useSecantHessVec_ ) {
302 if ( useSecantPrecond_ && !useSecantHessVec_ ) {
303 os <<
" with " <<
ESecantToString(esec_) <<
" Preconditioning" << std::endl;
305 else if ( !useSecantPrecond_ && useSecantHessVec_ ) {
306 os <<
" with " <<
ESecantToString(esec_) <<
" Hessian Approximation" << std::endl;
309 os <<
" with " <<
ESecantToString(esec_) <<
" Preconditioning and Hessian Approximation" << std::endl;
318template<
typename Real>
319void TrustRegionAlgorithm<Real>::writeOutput(std::ostream& os,
bool print_header)
const {
320 std::ios_base::fmtflags osFlags(os.flags());
321 os << std::scientific << std::setprecision(6);
322 if ( state_->iter == 0 ) {
325 if ( print_header ) {
328 if ( state_->iter == 0 ) {
330 os << std::setw(6) << std::left << state_->iter;
331 os << std::setw(15) << std::left << state_->value;
332 os << std::setw(15) << std::left << state_->gnorm;
333 os << std::setw(15) << std::left <<
"---";
334 os << std::setw(15) << std::left << state_->searchSize;
335 os << std::setw(10) << std::left << state_->nfval;
336 os << std::setw(10) << std::left << state_->ngrad;
337 os << std::setw(10) << std::left <<
"---";
338 if ( etr_ == TRUSTREGION_U_TRUNCATEDCG || etr_ == TRUSTREGION_U_SPG ) {
339 os << std::setw(10) << std::left <<
"---";
340 os << std::setw(10) << std::left <<
"---";
346 os << std::setw(6) << std::left << state_->iter;
347 os << std::setw(15) << std::left << state_->value;
348 os << std::setw(15) << std::left << state_->gnorm;
349 os << std::setw(15) << std::left << state_->snorm;
350 os << std::setw(15) << std::left << state_->searchSize;
351 os << std::setw(10) << std::left << state_->nfval;
352 os << std::setw(10) << std::left << state_->ngrad;
353 os << std::setw(10) << std::left << TRflag_;
354 if ( etr_ == TRUSTREGION_U_TRUNCATEDCG || etr_ == TRUSTREGION_U_SPG ) {
355 os << std::setw(10) << std::left << SPiter_;
356 os << 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)