ROL
ROL_TypeB_NewtonKrylovAlgorithm_Def.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_TYPEB_NEWTONKRYLOVALGORITHM_DEF_HPP
11#define ROL_TYPEB_NEWTONKRYLOVALGORITHM_DEF_HPP
12
13namespace ROL {
14namespace TypeB {
15
16template<typename Real>
18 const Ptr<Secant<Real>> &secant)
19 : secant_(secant), esec_(SECANT_USERDEFINED) {
21
22 if ( secant_ == nullPtr ) {
23 secantName_ = list.sublist("General").sublist("Secant").get("Type","Limited-Memory BFGS");
26 }
27 else {
28 secantName_ = list.sublist("General").sublist("Secant").get("User Defined Secant Name",
29 "Unspecified User Defined Secant Method");
30 }
31
32 krylovName_ = list.sublist("General").sublist("Krylov").get("Type","Conjugate Gradients");
35}
36
37template<typename Real>
39 const Ptr<Krylov<Real>> &krylov,
40 const Ptr<Secant<Real>> &secant)
43
44 if ( secant_ == nullPtr ) {
45 secantName_ = list.sublist("General").sublist("Secant").get("Type","Limited-Memory BFGS");
48 }
49 else {
50 secantName_ = list.sublist("General").sublist("Secant").get("User Defined Secant Name",
51 "Unspecified User Defined Secant Method");
52 }
53
54 krylovName_ = list.sublist("General").sublist("Krylov").get("User Defined Krylov Name",
55 "Unspecified User Defined Krylov Method");
56}
57
58template<typename Real>
60 // Set status test
61 status_->reset();
62 status_->add(makePtr<StatusTest<Real>>(list));
63
64 // Parse parameter list
65 ParameterList &lslist = list.sublist("Step").sublist("Line Search");
66 maxit_ = lslist.get("Function Evaluation Limit", 20);
67 alpha0_ = lslist.get("Initial Step Size", 1.0);
68 useralpha_ = lslist.get("User Defined Initial Step Size", false);
69 usePrevAlpha_ = lslist.get("Use Previous Step Length as Initial Guess", false);
70 c1_ = lslist.get("Sufficient Decrease Tolerance", 1e-4);
71 rhodec_ = lslist.sublist("Line-Search Method").get("Backtracking Rate", 0.5);
72
73 useSecantHessVec_ = list.sublist("General").sublist("Secant").get("Use as Hessian", false);
74 useSecantPrecond_ = list.sublist("General").sublist("Secant").get("Use as Preconditioner", false);
75
76 verbosity_ = list.sublist("General").get("Output Level", 0);
78}
79
80template<typename Real>
82 const Vector<Real> &g,
83 Objective<Real> &obj,
85 std::ostream &outStream) {
86 const Real one(1);
87 if (proj_ == nullPtr) {
88 proj_ = makePtr<PolyhedralProjection<Real>>(makePtrFromRef(bnd));
89 }
90 // Initialize data
92 // Update approximate gradient and approximate objective function.
93 Real ftol = std::sqrt(ROL_EPSILON<Real>());
94 proj_->project(x,outStream);
95 state_->iterateVec->set(x);
96 obj.update(x,UpdateType::Initial,state_->iter);
97 state_->value = obj.value(x,ftol); state_->nfval++;
98 obj.gradient(*state_->gradientVec,x,ftol); state_->ngrad++;
99 state_->stepVec->set(x);
100 state_->stepVec->axpy(-one,state_->gradientVec->dual());
101 proj_->project(*state_->stepVec,outStream);
102 state_->stepVec->axpy(-one,x);
103 state_->gnorm = state_->stepVec->norm();
104 state_->snorm = ROL_INF<Real>();
105 if (!useralpha_) alpha0_ = one;
106 state_->searchSize = alpha0_;
107}
108
109template<typename Real>
111 const Vector<Real> &g,
112 Objective<Real> &obj,
114 std::ostream &outStream ) {
115 const Real one(1);
116 // Initialize trust-region data
117 initialize(x,g,obj,bnd,outStream);
118 Ptr<Vector<Real>> s = x.clone(), gp = x.clone(), gold = g.clone();
119 Ptr<Vector<Real>> pwa = x.clone(), pwa1 = x.clone();
120 Real ftrial(0), gs(0), tol(std::sqrt(ROL_EPSILON<Real>()));
121
122 Ptr<LinearOperator<Real>> hessian, precond;
123
124 // Output
125 if (verbosity_ > 0) writeOutput(outStream,true);
126
127 // Compute steepest descent step
128 gp->set(state_->gradientVec->dual());
129 while (status_->check(*state_)) {
130 // Compute step
131 hessian = makePtr<HessianPNK>(makePtrFromRef(obj),makePtrFromRef(bnd),
132 state_->iterateVec,state_->gradientVec,state_->gnorm,
134 precond = makePtr<PrecondPNK>(makePtrFromRef(obj),makePtrFromRef(bnd),
135 state_->iterateVec,state_->gradientVec,state_->gnorm,
137 flagKrylov_ = 0;
138 krylov_->run(*s,*hessian,*state_->gradientVec,*precond,iterKrylov_,flagKrylov_);
139 if (flagKrylov_ == 2 && iterKrylov_ <= 1) {
140 s->set(*gp);
141 }
142 // Perform backtracking line search
143 if (!usePrevAlpha_) state_->searchSize = alpha0_;
144 x.set(*state_->iterateVec);
145 x.axpy(-state_->searchSize,*s);
146 proj_->project(x,outStream);
148 ftrial = obj.value(x,tol); ls_nfval_ = 1;
149 state_->stepVec->set(x);
150 state_->stepVec->axpy(-one,*state_->iterateVec);
151 gs = state_->stepVec->dot(*gp);
152 if (verbosity_ > 1) {
153 outStream << " In TypeB::NewtonKrylovAlgorithm: Line Search" << std::endl;
154 outStream << " Step size: " << state_->searchSize << std::endl;
155 outStream << " Trial objective value: " << ftrial << std::endl;
156 outStream << " Computed reduction: " << state_->value-ftrial << std::endl;
157 outStream << " Dot product of gradient and step: " << gs << std::endl;
158 outStream << " Sufficient decrease bound: " << -gs*c1_ << std::endl;
159 outStream << " Number of function evaluations: " << ls_nfval_ << std::endl;
160 }
161 while ( state_->value - ftrial < -c1_*gs && ls_nfval_ < maxit_ ) {
162 state_->searchSize *= rhodec_;
163 x.set(*state_->iterateVec);
164 x.axpy(-state_->searchSize,*s);
165 proj_->project(x,outStream);
167 ftrial = obj.value(x,tol); ls_nfval_++;
168 state_->stepVec->set(x);
169 state_->stepVec->axpy(-one,*state_->iterateVec);
170 gs = state_->stepVec->dot(*gp);
171 if (verbosity_ > 1) {
172 outStream << std::endl;
173 outStream << " Step size: " << state_->searchSize << std::endl;
174 outStream << " Trial objective value: " << ftrial << std::endl;
175 outStream << " Computed reduction: " << state_->value-ftrial << std::endl;
176 outStream << " Dot product of gradient and step: " << gs << std::endl;
177 outStream << " Sufficient decrease bound: " << -gs*c1_ << std::endl;
178 outStream << " Number of function evaluations: " << ls_nfval_ << std::endl;
179 }
180 }
181 state_->nfval += ls_nfval_;
182
183 // Compute norm of step
184 state_->snorm = state_->stepVec->norm();
185
186 // Update iterate
187 state_->iterateVec->set(x);
188
189 // Compute new value and gradient
190 state_->iter++;
191 state_->value = ftrial;
192 obj.update(x,UpdateType::Accept,state_->iter);
193 gold->set(*state_->gradientVec);
194 obj.gradient(*state_->gradientVec,x,tol); state_->ngrad++;
195 gp->set(state_->gradientVec->dual());
196
197 // Compute projected gradient norm
198 s->set(x); s->axpy(-one,*gp);
199 proj_->project(*s,outStream);
200 s->axpy(-one,x);
201 state_->gnorm = s->norm();
202
203 // Update secant
204 secant_->updateStorage(x,*state_->gradientVec,*gold,*state_->stepVec,state_->snorm,state_->iter);
205
206 // Update Output
207 if (verbosity_ > 0) writeOutput(outStream,writeHeader_);
208 }
210}
211
212template<typename Real>
214 std::ostream &outStream ) {
215 if (problem.getPolyhedralProjection() == nullPtr) {
216 return TypeB::Algorithm<Real>::run(problem,outStream);
217 }
218 else {
219 throw Exception::NotImplemented(">>> TypeB::NewtonKrylovAlgorithm::run : This algorithm cannot solve problems with linear equality constraints!");
220 }
221}
222
223template<typename Real>
225 const Vector<Real> &g,
226 Objective<Real> &obj,
228 Constraint<Real> &linear_econ,
229 Vector<Real> &linear_emul,
230 const Vector<Real> &linear_eres,
231 std::ostream &outStream ) {
232 throw Exception::NotImplemented(">>> TypeB::NewtonKrylovAlgorithm::run : This algorithm cannot solve problems with linear equality constraints!");
233}
234
235template<typename Real>
236void NewtonKrylovAlgorithm<Real>::writeHeader( std::ostream& os ) const {
237 std::ios_base::fmtflags osFlags(os.flags());
238 if (verbosity_ > 1) {
239 os << std::string(114,'-') << std::endl;
240 if (!useSecantHessVec_) {
241 os << "Line-Search Projected Newton";
242 }
243 else {
244 os << "Line-Search Projected Quasi-Newton with " << secantName_ << " Hessian approximation";
245 }
246 os << " 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 << " alpha - Line search step length" << std::endl;
252 os << " #fval - Cumulative number of times the objective function was evaluated" << std::endl;
253 os << " #grad - Cumulative number of times the gradient was computed" << std::endl;
254 os << " ls_#fval - Number of the times the objective function was evaluated during the line search" << std::endl;
255 os << " iterCG - Number of Krylov iterations" << std::endl << std::endl;
256 os << " flagGC - Krylov flag" << std::endl;
257 for( int flag = CG_FLAG_SUCCESS; flag != CG_FLAG_UNDEFINED; ++flag ) {
258 os << " " << NumberToString(flag) << " - "
259 << ECGFlagToString(static_cast<ECGFlag>(flag)) << std::endl;
260 }
261 os << std::string(114,'-') << std::endl;
262 }
263
264 os << " ";
265 os << std::setw(6) << std::left << "iter";
266 os << std::setw(15) << std::left << "value";
267 os << std::setw(15) << std::left << "gnorm";
268 os << std::setw(15) << std::left << "snorm";
269 os << std::setw(15) << std::left << "alpha";
270 os << std::setw(10) << std::left << "#fval";
271 os << std::setw(10) << std::left << "#grad";
272 os << std::setw(10) << std::left << "#ls_fval";
273 os << std::setw(10) << std::left << "iterCG";
274 os << std::setw(10) << std::left << "flagCG";
275 os << std::endl;
276 os.flags(osFlags);
277}
278
279template<typename Real>
280void NewtonKrylovAlgorithm<Real>::writeName( std::ostream& os ) const {
281 std::ios_base::fmtflags osFlags(os.flags());
282 if (!useSecantHessVec_) {
283 os << std::endl << "Line-Search Projected Newton (Type B, Bound Constraints)" << std::endl;
284 }
285 else {
286 os << std::endl << "Line-Search Projected Quasi-Newton with "
287 << secantName_ << " Hessian approximation" << std::endl;
288 }
289 os.flags(osFlags);
290}
291
292template<typename Real>
293void NewtonKrylovAlgorithm<Real>::writeOutput( std::ostream& os, bool write_header ) const {
294 std::ios_base::fmtflags osFlags(os.flags());
295 os << std::scientific << std::setprecision(6);
296 if ( state_->iter == 0 ) writeName(os);
297 if ( write_header ) writeHeader(os);
298 if ( state_->iter == 0 ) {
299 os << " ";
300 os << std::setw(6) << std::left << state_->iter;
301 os << std::setw(15) << std::left << state_->value;
302 os << std::setw(15) << std::left << state_->gnorm;
303 os << std::setw(15) << std::left << "---";
304 os << std::setw(15) << std::left << "---";
305 os << std::setw(10) << std::left << state_->nfval;
306 os << std::setw(10) << std::left << state_->ngrad;
307 os << std::setw(10) << std::left << "---";
308 os << std::setw(10) << std::left << "---";
309 os << std::setw(10) << std::left << "---";
310 os << std::endl;
311 }
312 else {
313 os << " ";
314 os << std::setw(6) << std::left << state_->iter;
315 os << std::setw(15) << std::left << state_->value;
316 os << std::setw(15) << std::left << state_->gnorm;
317 os << std::setw(15) << std::left << state_->snorm;
318 os << std::setw(15) << std::left << state_->searchSize;
319 os << std::setw(10) << std::left << state_->nfval;
320 os << std::setw(10) << std::left << state_->ngrad;
321 os << std::setw(10) << std::left << ls_nfval_;
322 os << std::setw(10) << std::left << iterKrylov_;
323 os << std::setw(10) << std::left << flagKrylov_;
324 os << std::endl;
325 }
326 os.flags(osFlags);
327}
328
329} // namespace TypeB
330} // namespace ROL
331
332#endif
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.
Provides the interface to apply upper and lower bound constraints.
Defines the general constraint operator interface.
Provides definitions for Krylov solvers.
Provides the interface to evaluate objective functions.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
const Ptr< PolyhedralProjection< Real > > & getPolyhedralProjection()
Get the polyhedral projection object. This is a null pointer if no linear constraints and/or bounds a...
Provides interface for and implements limited-memory secant operators.
Provides an interface to check status of optimization algorithms.
virtual void run(Problem< Real > &problem, std::ostream &outStream=std::cout)
Run algorithm on bound constrained problems (Type-B). This is the primary Type-B interface.
Ptr< PolyhedralProjection< Real > > proj_
void initialize(const Vector< Real > &x, const Vector< Real > &g)
virtual void writeExitStatus(std::ostream &os) const
const Ptr< AlgorithmState< Real > > state_
const Ptr< CombinedStatusTest< Real > > status_
Ptr< Secant< Real > > secant_
Secant object (used for quasi-Newton).
Real alpha0_
Initial line search step size (default: 1.0).
void run(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout) override
Run algorithm on bound constrained problems (Type-B). This general interface supports the use of dual...
bool useSecantPrecond_
Whether or not to use a secant approximation to precondition inexact Newton.
int maxit_
Maximum number of line search steps (default: 20).
int flagKrylov_
Termination flag for Krylov method (used for inexact Newton).
bool useralpha_
Flag to use user-defined initial step size (default: false).
NewtonKrylovAlgorithm(ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr)
Real c1_
Sufficient Decrease Parameter (default: 1e-4).
bool useSecantHessVec_
Whether or not to use to a secant approximation as the Hessian.
void writeName(std::ostream &os) const override
Print step name.
Real rhodec_
Backtracking rate (default: 0.5).
int iterKrylov_
Number of Krylov iterations (used for inexact Newton).
Ptr< Krylov< Real > > krylov_
Krylov solver object (used for inexact Newton).
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout)
bool usePrevAlpha_
Flag to use previous step size as next initial step size (default: false).
void writeHeader(std::ostream &os) const override
Print iterate header.
void writeOutput(std::ostream &os, const bool write_header=false) const override
Print iterate status.
Defines the linear algebra or vector space interface.
virtual void set(const Vector &x)
Set where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
EKrylov StringToEKrylov(std::string s)
std::string NumberToString(T Number)
Definition ROL_Types.hpp:47
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
Definition ROL_Types.hpp:57
ESecant StringToESecant(std::string s)
Ptr< Krylov< Real > > KrylovFactory(ParameterList &parlist)
@ SECANT_USERDEFINED
@ CG_FLAG_UNDEFINED
@ CG_FLAG_SUCCESS
ROL::Ptr< Secant< Real > > SecantFactory(ROL::ParameterList &parlist, ESecantMode mode=SECANTMODE_BOTH)
std::string ECGFlagToString(ECGFlag cgf)
Real ROL_INF(void)
Definition ROL_Types.hpp:71