ROL
ROL_TypeP_TrustRegionAlgorithm.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_TYPEP_TRUSTREGIONALGORITHM_HPP
11#define ROL_TYPEP_TRUSTREGIONALGORITHM_HPP
12
16#include "ROL_Types.hpp"
17
22namespace ROL {
23namespace TypeP {
24
31
32inline std::string ETrustRegionPToString(ETrustRegionP alg) {
33 std::string retString;
34 switch(alg) {
35 case TRUSTREGION_P_SPG: retString = "SPG"; break;
36 case TRUSTREGION_P_SPG2: retString = "Simplified SPG"; break;
37 case TRUSTREGION_P_NCG: retString = "NCG"; break;
38 case TRUSTREGION_P_LAST: retString = "Last Type (Dummy)"; break;
39 default: retString = "INVALID ETrustRegionP";
40 }
41 return retString;
42}
43
45 return( (alg == TRUSTREGION_P_SPG) ||
46 (alg == TRUSTREGION_P_SPG2) ||
47 (alg == TRUSTREGION_P_NCG) ||
48 (alg == TRUSTREGION_P_LAST)
49 );
50}
51
53 return type = static_cast<ETrustRegionP>(type+1);
54}
55
57 ETrustRegionP oldval = type;
58 ++type;
59 return oldval;
60}
61
63 return type = static_cast<ETrustRegionP>(type-1);
64}
65
67 ETrustRegionP oldval = type;
68 --type;
69 return oldval;
70}
71
72inline ETrustRegionP StringToETrustRegionP(std::string s) {
73 s = removeStringFormat(s);
74 for ( ETrustRegionP alg = TRUSTREGION_P_SPG; alg < TRUSTREGION_P_LAST; alg++ ) {
75 if ( !s.compare(removeStringFormat(ETrustRegionPToString(alg))) ) {
76 return alg;
77 }
78 }
79 return TRUSTREGION_P_SPG;
80}
81
82template<typename Real>
84private:
85 Ptr<TrustRegionModel_U<Real>> model_;
86
87 // TRUST REGION PARAMETERS
88 Real delMax_;
89 Real eta0_;
90 Real eta1_;
91 Real eta2_;
92 Real gamma0_;
93 Real gamma1_;
94 Real gamma2_;
95 Real TRsafe_;
96 Real eps_;
98
99 // ITERATION FLAGS/INFORMATION
103
104 // SECANT INFORMATION
108
109 // TRUNCATED CG INFORMATION
110 Real tol1_;
111 Real tol2_;
112 int maxit_;
113
114 // ALGORITHM SPECIFIC PARAMETERS
115 bool useNM_;
117 Real mu0_;
118 Real spexp_;
121 Real alpha_;
123 Real interpf_;
124 Real extrapf_;
125 Real qtol_;
128 Real gamma_;
134 Real etaNCG_;
136
137 // Inexactness Parameters
138 std::vector<bool> useInexact_;
141 Real scale_;
142 Real omega_;
143 Real force_;
146 Real gtol_;
147
149 Real t0_;
150
151 mutable int nhess_;
152 unsigned verbosity_;
154
155 using TypeP::Algorithm<Real>::state_;
156 using TypeP::Algorithm<Real>::status_;
157 using TypeP::Algorithm<Real>::pgstep;
158
159 void initialize(Vector<Real> &x,
160 const Vector<Real> &g,
161 Real ftol,
162 Objective<Real> &sobj,
163 Objective<Real> &nobj,
164 Vector<Real> &px,
165 Vector<Real> &dg,
166 std::ostream &outStream = std::cout);
167
168public:
169 TrustRegionAlgorithm(ParameterList &list, const Ptr<Secant<Real>> &secant = nullPtr);
170
171 using TypeP::Algorithm<Real>::run;
172 void run( Vector<Real> &x,
173 const Vector<Real> &g,
174 Objective<Real> &sobj,
175 Objective<Real> &nobj,
176 std::ostream &outStream = std::cout) override;
177
178 void writeHeader( std::ostream& os ) const override;
179
180 void writeName( std::ostream& os ) const override;
181
182 void writeOutput( std::ostream& os, bool write_header = false ) const override;
183
184private:
185
186 Real computeValue(Real inTol,
187 Real &outTol,
188 Real pRed,
189 Real &fold,
190 int iter,
191 const Vector<Real> &x,
192 const Vector<Real> &xold,
193 Objective<Real> &obj);
194
195 void computeGradient(const Vector<Real> &x,
196 Vector<Real> &g,
197 Vector<Real> &px,
198 Vector<Real> &dg,
199 Vector<Real> &pwa,
200 Real del,
201 Objective<Real> &sobj,
202 Objective<Real> &nobj,
203 bool accept,
204 Real &gtol,
205 Real &gnorm,
206 std::ostream &outStream = std::cout) const;
207
208 // Compute the projected step s = P(x + alpha*w) - x
209 // Returns the norm of the projected step s
210 // s -- The projected step upon return
211 // w -- The direction vector w (unchanged)
212 // x -- The anchor vector x (unchanged)
213 // alpha -- The step size (unchanged)
214 // Real dgpstep(Vector<Real> &s, const Vector<Real> &w,
215 // const Vector<Real> &x, const Real alpha,
216 // std::ostream &outStream = std::cout) const;
217
218 // Compute Cauchy point, i.e., the minimizer of q(P(x - alpha*g)-x)
219 // subject to the trust region constraint ||P(x - alpha*g)-x|| <= del
220 // s -- The Cauchy step upon return: Primal optimization space vector
221 // alpha -- The step length for the Cauchy point upon return
222 // x -- The anchor vector x (unchanged): Primal optimization space vector
223 // g -- The (dual) gradient vector g (unchanged): Primal optimization space vector
224 // del -- The trust region radius (unchanged)
225 // model -- Trust region model
226 // dwa -- Dual working array, stores Hessian applied to step
227 // dwa1 -- Dual working array
228 Real dcauchy(Vector<Real> &s,
229 Real &alpha,
230 Real &sval,
231 Real &nval,
232 const Vector<Real> &x,
233 const Vector<Real> &g,
234 Real del,
236 Objective<Real> &nobj,
237 Vector<Real> &px,
238 Vector<Real> &dwa,
239 Vector<Real> &dwa1,
240 std::ostream &outStream = std::cout);
241
242 void dspg2(Vector<Real> &y,
243 Real &sval,
244 Real &nval,
245 Real &pRed,
246 Vector<Real> &gmod,
247 const Vector<Real> &x,
248 Real del,
250 Objective<Real> &nobj,
251 Vector<Real> &pwa,
252 Vector<Real> &pwa1,
253 Vector<Real> &pwa2,
254 Vector<Real> &dwa,
255 std::ostream &outStream = std::cout);
256
257 void dspg(Vector<Real> &y,
258 Real &sval,
259 Real &nval,
260 Vector<Real> &gmod,
261 const Vector<Real> &x,
262 Real del,
264 Objective<Real> &nobj,
265 Vector<Real> &ymin,
266 Vector<Real> &pwa,
267 Vector<Real> &pwa1,
268 Vector<Real> &pwa2,
269 Vector<Real> &pwa3,
270 Vector<Real> &pwa4,
271 Vector<Real> &pwa5,
272 Vector<Real> &dwa,
273 std::ostream &outStream = std::cout);
274
275 void dncg(Vector<Real> &y,
276 Real &sval,
277 Real &nval,
278 Vector<Real> &gmod,
279 const Vector<Real> &x,
280 Real del,
282 Objective<Real> &nobj,
283 Vector<Real> &s,
284 Vector<Real> &pwa1,
285 Vector<Real> &pwa2,
286 Vector<Real> &pwa3,
287 Vector<Real> &pwa4,
288 Vector<Real> &pwa5,
289 Vector<Real> &dwa,
290 std::ostream &outStream = std::cout);
291
292
293 void dprox(Vector<Real> &x,
294 const Vector<Real> &x0,
295 Real t,
296 Real del,
297 Objective<Real> &nobj,
298 Vector<Real> &y0,
299 Vector<Real> &y1,
300 Vector<Real> &yc,
301 Vector<Real> &pwa,
302 std::ostream &outStream = std::cout) const;
303
304 void dbls(Real &alpha, Real &nval, Real &pred,
305 const Vector<Real> &y,
306 const Vector<Real> &s,
307 Real lambda, Real tmax,
308 Real kappa, Real gs,
309 Objective<Real> &nobj,
310 Vector<Real> &pwa);
311
312}; // class ROL::TypeP::TrustRegionAlgorithm
313
314} // namespace TypeP
315} // namespace ROL
316
318
319#endif
Contains definitions of custom data types in ROL.
Provides the interface to evaluate objective functions.
Provides interface for and implements limited-memory secant operators.
Provides the interface to evaluate trust-region model functions.
Provides an interface to run optimization algorithms to minimize composite optimization problems f+ph...
void pgstep(Vector< Real > &pgiter, Vector< Real > &pgstep, Objective< Real > &nobj, const Vector< Real > &x, const Vector< Real > &dg, Real t, Real &tol) const
const Ptr< AlgorithmState< Real > > state_
const Ptr< CombinedStatusTest< Real > > status_
Provides an interface to run the trust-region algorithm.
TRUtils::ETRFlag TRflag_
Trust-region exit flag.
void dncg(Vector< Real > &y, Real &sval, Real &nval, Vector< Real > &gmod, const Vector< Real > &x, Real del, TrustRegionModel_U< Real > &model, Objective< Real > &nobj, Vector< Real > &s, Vector< Real > &pwa1, Vector< Real > &pwa2, Vector< Real > &pwa3, Vector< Real > &pwa4, Vector< Real > &pwa5, Vector< Real > &dwa, std::ostream &outStream=std::cout)
unsigned verbosity_
Output level (default: 0)
ESecant esec_
Secant type (default: Limited-Memory BFGS)
Real gamma1_
Radius decrease rate (positive rho) (default: 0.25)
void dprox(Vector< Real > &x, const Vector< Real > &x0, Real t, Real del, Objective< Real > &nobj, Vector< Real > &y0, Vector< Real > &y1, Vector< Real > &yc, Vector< Real > &pwa, std::ostream &outStream=std::cout) const
void dspg(Vector< Real > &y, Real &sval, Real &nval, Vector< Real > &gmod, const Vector< Real > &x, Real del, TrustRegionModel_U< Real > &model, Objective< Real > &nobj, Vector< Real > &ymin, Vector< Real > &pwa, Vector< Real > &pwa1, Vector< Real > &pwa2, Vector< Real > &pwa3, Vector< Real > &pwa4, Vector< Real > &pwa5, Vector< Real > &dwa, std::ostream &outStream=std::cout)
void initialize(Vector< Real > &x, const Vector< Real > &g, Real ftol, Objective< Real > &sobj, Objective< Real > &nobj, Vector< Real > &px, Vector< Real > &dg, std::ostream &outStream=std::cout)
void writeName(std::ostream &os) const override
Print step name.
Real gamma0_
Radius decrease rate (negative rho) (default: 0.0625)
void writeOutput(std::ostream &os, bool write_header=false) const override
Print iterate status.
bool interpRad_
Interpolate the trust-region radius if ratio is negative (default: false)
int SPiter_
Subproblem solver iteration count.
int explim_
Maximum number of Cauchy point expansion steps (default: 10)
Real computeValue(Real inTol, Real &outTol, Real pRed, Real &fold, int iter, const Vector< Real > &x, const Vector< Real > &xold, Objective< Real > &obj)
Real tol1_
Absolute tolerance for truncated CG (default: 1e-4)
Real mu0_
Sufficient decrease parameter (default: 1e-2)
Real qtol_
Relative tolerance for computed decrease in Cauchy point computation (default: 1-8)
void dbls(Real &alpha, Real &nval, Real &pred, const Vector< Real > &y, const Vector< Real > &s, Real lambda, Real tmax, Real kappa, Real gs, Objective< Real > &nobj, Vector< Real > &pwa)
Real interpf_
Backtracking rate for Cauchy point computation (default: 1e-1)
Real delMax_
Maximum trust-region radius (default: ROL_INF)
Real alpha_
Initial Cauchy point step length (default: 1.0)
Real dcauchy(Vector< Real > &s, Real &alpha, Real &sval, Real &nval, const Vector< Real > &x, const Vector< Real > &g, Real del, TrustRegionModel_U< Real > &model, Objective< Real > &nobj, Vector< Real > &px, Vector< Real > &dwa, Vector< Real > &dwa1, std::ostream &outStream=std::cout)
bool normAlpha_
Normalize initial Cauchy point step length (default: false)
Real TRsafe_
Safeguard size for numerically evaluating ratio (default: 1e2)
void computeGradient(const Vector< Real > &x, Vector< Real > &g, Vector< Real > &px, Vector< Real > &dg, Vector< Real > &pwa, Real del, Objective< Real > &sobj, Objective< Real > &nobj, bool accept, Real &gtol, Real &gnorm, std::ostream &outStream=std::cout) const
Real extrapf_
Extrapolation rate for Cauchy point computation (default: 1e1)
Real eta2_
Radius increase threshold (default: 0.9)
void writeHeader(std::ostream &os) const override
Print iterate header.
Ptr< TrustRegionModel_U< Real > > model_
Container for trust-region model.
Real eps_
Safeguard for numerically evaluating ratio.
int redlim_
Maximum number of Cauchy point reduction steps (default: 10)
void run(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &sobj, Objective< Real > &nobj, std::ostream &outStream=std::cout) override
Run algorithm on unconstrained problems (Type-U). This general interface supports the use of dual opt...
bool writeHeader_
Flag to write header at every iteration.
Real tol2_
Relative tolerance for truncated CG (default: 1e-2)
Real spexp_
Relative tolerance exponent for subproblem solve (default: 1, range: [1,2])
Real eta1_
Radius decrease threshold (default: 0.05)
int SPflag_
Subproblem solver termination flag.
bool useSecantHessVec_
Flag to use secant as Hessian (default: false)
Real eta0_
Step acceptance threshold (default: 0.05)
int maxit_
Maximum number of CG iterations (default: 25)
bool useSecantPrecond_
Flag to use secant as a preconditioner (default: false)
Real gamma2_
Radius increase rate (default: 2.5)
int nhess_
Number of Hessian applications.
void dspg2(Vector< Real > &y, Real &sval, Real &nval, Real &pRed, Vector< Real > &gmod, const Vector< Real > &x, Real del, TrustRegionModel_U< Real > &model, Objective< Real > &nobj, Vector< Real > &pwa, Vector< Real > &pwa1, Vector< Real > &pwa2, Vector< Real > &dwa, std::ostream &outStream=std::cout)
Defines the linear algebra or vector space interface.
ETrustRegionP StringToETrustRegionP(std::string s)
std::string ETrustRegionPToString(ETrustRegionP alg)
EAlgorithmP & operator--(EAlgorithmP &type)
EAlgorithmP & operator++(EAlgorithmP &type)
int isValidTrustRegionP(ETrustRegionP alg)
std::string removeStringFormat(std::string s)