ROL
ROL_ProjectedNewtonKrylovStep.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_PROJECTEDNEWTONKRYLOVSTEP_H
11#define ROL_PROJECTEDNEWTONKRYLOVSTEP_H
12
13#include "ROL_Types.hpp"
14#include "ROL_Step.hpp"
15
16#include "ROL_Secant.hpp"
17#include "ROL_KrylovFactory.hpp"
19
20#include <sstream>
21#include <iomanip>
22
30namespace ROL {
31
32template <class Real>
33class ProjectedNewtonKrylovStep : public Step<Real> {
34private:
35
36 ROL::Ptr<Secant<Real> > secant_;
37 ROL::Ptr<Krylov<Real> > krylov_;
38
41
42 ROL::Ptr<Vector<Real> > gp_;
43 ROL::Ptr<Vector<Real> > d_;
44
48 const bool computeObj_;
49
52
53 std::string krylovName_;
54 std::string secantName_;
55
56 class HessianPNK : public LinearOperator<Real> {
57 private:
58 const ROL::Ptr<Objective<Real> > obj_;
59 const ROL::Ptr<BoundConstraint<Real> > bnd_;
60 const ROL::Ptr<Vector<Real> > x_;
61 const ROL::Ptr<Vector<Real> > g_;
62 ROL::Ptr<Vector<Real> > v_;
63 Real eps_;
64 public:
65 HessianPNK(const ROL::Ptr<Objective<Real> > &obj,
66 const ROL::Ptr<BoundConstraint<Real> > &bnd,
67 const ROL::Ptr<Vector<Real> > &x,
68 const ROL::Ptr<Vector<Real> > &g,
69 Real eps = 0 )
70 : obj_(obj), bnd_(bnd), x_(x), g_(g), eps_(eps) {
71 v_ = x_->clone();
72 }
73 void apply(Vector<Real> &Hv, const Vector<Real> &v, Real &tol) const {
74 v_->set(v);
75 bnd_->pruneActive(*v_,*g_,*x_,eps_);
76 obj_->hessVec(Hv,*v_,*x_,tol);
77 bnd_->pruneActive(Hv,*g_,*x_,eps_);
78 v_->set(v);
79 bnd_->pruneInactive(*v_,*g_,*x_,eps_);
80 Hv.plus(v_->dual());
81 }
82 };
83
84 class PrecondPNK : public LinearOperator<Real> {
85 private:
86 const ROL::Ptr<Objective<Real> > obj_;
87 const ROL::Ptr<Secant<Real> > secant_;
88 const ROL::Ptr<BoundConstraint<Real> > bnd_;
89 const ROL::Ptr<Vector<Real> > x_;
90 const ROL::Ptr<Vector<Real> > g_;
91 ROL::Ptr<Vector<Real> > v_;
92 Real eps_;
93 const bool useSecant_;
94 public:
95 PrecondPNK(const ROL::Ptr<Objective<Real> > &obj,
96 const ROL::Ptr<BoundConstraint<Real> > &bnd,
97 const ROL::Ptr<Vector<Real> > &x,
98 const ROL::Ptr<Vector<Real> > &g,
99 Real eps = 0 )
100 : obj_(obj), bnd_(bnd), x_(x), g_(g), eps_(eps), useSecant_(false) {
101 v_ = x_->clone();
102 }
103 PrecondPNK(const ROL::Ptr<Secant<Real> > &secant,
104 const ROL::Ptr<BoundConstraint<Real> > &bnd,
105 const ROL::Ptr<Vector<Real> > &x,
106 const ROL::Ptr<Vector<Real> > &g,
107 Real eps = 0 )
108 : secant_(secant), bnd_(bnd), x_(x), g_(g), eps_(eps), useSecant_(true) {
109 v_ = x_->clone();
110 }
111 void apply(Vector<Real> &Hv, const Vector<Real> &v, Real &tol) const {
112 Hv.set(v.dual());
113 }
114 void applyInverse(Vector<Real> &Hv, const Vector<Real> &v, Real &tol) const {
115 v_->set(v);
116 bnd_->pruneActive(*v_,*g_,*x_,eps_);
117 if ( useSecant_ ) {
118 secant_->applyH(Hv,*v_);
119 }
120 else {
121 obj_->precond(Hv,*v_,*x_,tol);
122 }
123 bnd_->pruneActive(Hv,*g_,*x_,eps_);
124 v_->set(v);
125 bnd_->pruneInactive(*v_,*g_,*x_,eps_);
126 Hv.plus(v_->dual());
127 }
128 };
129
130public:
131
132 using Step<Real>::initialize;
133 using Step<Real>::compute;
134 using Step<Real>::update;
135
143 ProjectedNewtonKrylovStep( ROL::ParameterList &parlist, const bool computeObj = true )
144 : Step<Real>(), secant_(ROL::nullPtr), krylov_(ROL::nullPtr),
145 gp_(ROL::nullPtr), d_(ROL::nullPtr),
147 computeObj_(computeObj), useSecantPrecond_(false) {
148 // Parse ParameterList
149 ROL::ParameterList& Glist = parlist.sublist("General");
150 useSecantPrecond_ = Glist.sublist("Secant").get("Use as Preconditioner", false);
151 useProjectedGrad_ = Glist.get("Projected Gradient Criticality Measure", false);
152 verbosity_ = Glist.get("Print Verbosity",0);
153 // Initialize Krylov object
154 krylovName_ = Glist.sublist("Krylov").get("Type","Conjugate Gradients");
156 krylov_ = KrylovFactory<Real>(parlist);
157 // Initialize secant object
158 secantName_ = Glist.sublist("Secant").get("Type","Limited-Memory BFGS");
160 if ( useSecantPrecond_ ) {
161 secant_ = SecantFactory<Real>(parlist);
162 }
163 }
164
175 ProjectedNewtonKrylovStep(ROL::ParameterList &parlist,
176 const ROL::Ptr<Krylov<Real> > &krylov,
177 const ROL::Ptr<Secant<Real> > &secant,
178 const bool computeObj = true)
179 : Step<Real>(), secant_(secant), krylov_(krylov),
181 gp_(ROL::nullPtr), d_(ROL::nullPtr),
183 computeObj_(computeObj), useSecantPrecond_(false) {
184 // Parse ParameterList
185 ROL::ParameterList& Glist = parlist.sublist("General");
186 useSecantPrecond_ = Glist.sublist("Secant").get("Use as Preconditioner", false);
187 useProjectedGrad_ = Glist.get("Projected Gradient Criticality Measure", false);
188 verbosity_ = Glist.get("Print Verbosity",0);
189 // Initialize secant object
190 if ( useSecantPrecond_ ) {
191 if (secant_ == ROL::nullPtr ) {
192 secantName_ = Glist.sublist("Secant").get("Type","Limited-Memory BFGS");
194 secant_ = SecantFactory<Real>(parlist);
195 }
196 else {
197 secantName_ = Glist.sublist("Secant").get("User Defined Secant Name",
198 "Unspecified User Defined Secant Method");
199 }
200 }
201 // Initialize Krylov object
202 if ( krylov_ == ROL::nullPtr ) {
203 krylovName_ = Glist.sublist("Krylov").get("Type","Conjugate Gradients");
205 krylov_ = KrylovFactory<Real>(parlist);
206 }
207 }
208
209 void initialize( Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g,
211 AlgorithmState<Real> &algo_state ) {
212 Step<Real>::initialize(x,s,g,obj,bnd,algo_state);
213 gp_ = g.clone();
214 d_ = s.clone();
215 }
216
219 AlgorithmState<Real> &algo_state ) {
220 Real one(1);
221 ROL::Ptr<StepState<Real> > step_state = Step<Real>::getState();
222
223 // Build Hessian and Preconditioner object
224 ROL::Ptr<Objective<Real> > obj_ptr = ROL::makePtrFromRef(obj);
225 ROL::Ptr<BoundConstraint<Real> > bnd_ptr = ROL::makePtrFromRef(bnd);
226 ROL::Ptr<LinearOperator<Real> > hessian
227 = ROL::makePtr<HessianPNK>(obj_ptr,bnd_ptr,algo_state.iterateVec,
228 step_state->gradientVec,algo_state.gnorm);
229 ROL::Ptr<LinearOperator<Real> > precond;
230 if (useSecantPrecond_) {
231 precond = ROL::makePtr<PrecondPNK>(secant_,bnd_ptr,
232 algo_state.iterateVec,step_state->gradientVec,algo_state.gnorm);
233 }
234 else {
235 precond = ROL::makePtr<PrecondPNK>(obj_ptr,bnd_ptr,
236 algo_state.iterateVec,step_state->gradientVec,algo_state.gnorm);
237 }
238
239 // Run Krylov method
240 flagKrylov_ = 0;
241 krylov_->run(s,*hessian,*(step_state->gradientVec),*precond,iterKrylov_,flagKrylov_);
242
243 // Check Krylov flags
244 if ( flagKrylov_ == 2 && iterKrylov_ <= 1 ) {
245 s.set((step_state->gradientVec)->dual());
246 }
247 s.scale(-one);
248 }
249
250 void update( Vector<Real> &x, const Vector<Real> &s,
252 AlgorithmState<Real> &algo_state ) {
253 Real tol = std::sqrt(ROL_EPSILON<Real>()), one(1);
254 ROL::Ptr<StepState<Real> > step_state = Step<Real>::getState();
255 step_state->SPiter = iterKrylov_;
256 step_state->SPflag = flagKrylov_;
257
258 // Update iterate and store previous step
259 algo_state.iter++;
260 d_->set(x);
261 x.plus(s);
262 bnd.project(x);
263 (step_state->descentVec)->set(x);
264 (step_state->descentVec)->axpy(-one,*d_);
265 algo_state.snorm = s.norm();
266
267 // Compute new gradient
268 if ( useSecantPrecond_ ) {
269 gp_->set(*(step_state->gradientVec));
270 }
271 obj.update(x,true,algo_state.iter);
272 if ( computeObj_ ) {
273 algo_state.value = obj.value(x,tol);
274 algo_state.nfval++;
275 }
276 obj.gradient(*(step_state->gradientVec),x,tol);
277 algo_state.ngrad++;
278
279 // Update Secant Information
280 if ( useSecantPrecond_ ) {
281 secant_->updateStorage(x,*(step_state->gradientVec),*gp_,s,algo_state.snorm,algo_state.iter+1);
282 }
283
284 // Update algorithm state
285 (algo_state.iterateVec)->set(x);
286 if ( useProjectedGrad_ ) {
287 gp_->set(*(step_state->gradientVec));
288 bnd.computeProjectedGradient( *gp_, x );
289 algo_state.gnorm = gp_->norm();
290 }
291 else {
292 d_->set(x);
293 d_->axpy(-one,(step_state->gradientVec)->dual());
294 bnd.project(*d_);
295 d_->axpy(-one,x);
296 algo_state.gnorm = d_->norm();
297 }
298 }
299
300 std::string printHeader( void ) const {
301 std::stringstream hist;
302
303 if( verbosity_>0 ) {
304 hist << std::string(109,'-') << "\n";
306 hist << " status output definitions\n\n";
307 hist << " iter - Number of iterates (steps taken) \n";
308 hist << " value - Objective function value \n";
309 hist << " gnorm - Norm of the gradient\n";
310 hist << " snorm - Norm of the step (update to optimization vector)\n";
311 hist << " #fval - Cumulative number of times the objective function was evaluated\n";
312 hist << " #grad - Number of times the gradient was computed\n";
313 hist << " iterCG - Number of Krylov iterations used to compute search direction\n";
314 hist << " flagCG - Krylov solver flag" << "\n";
315 hist << std::string(109,'-') << "\n";
316 }
317
318 hist << " ";
319 hist << std::setw(6) << std::left << "iter";
320 hist << std::setw(15) << std::left << "value";
321 hist << std::setw(15) << std::left << "gnorm";
322 hist << std::setw(15) << std::left << "snorm";
323 hist << std::setw(10) << std::left << "#fval";
324 hist << std::setw(10) << std::left << "#grad";
325 hist << std::setw(10) << std::left << "iterCG";
326 hist << std::setw(10) << std::left << "flagCG";
327 hist << "\n";
328 return hist.str();
329 }
330 std::string printName( void ) const {
331 std::stringstream hist;
332 hist << "\n" << EDescentToString(DESCENT_NEWTONKRYLOV);
333 hist << " using " << krylovName_;
334 if ( useSecantPrecond_ ) {
335 hist << " with " << secantName_ << " preconditioning";
336 }
337 hist << "\n";
338 return hist.str();
339 }
340 std::string print( AlgorithmState<Real> &algo_state, bool print_header = false ) const {
341 std::stringstream hist;
342 hist << std::scientific << std::setprecision(6);
343 if ( algo_state.iter == 0 ) {
344 hist << printName();
345 }
346 if ( print_header ) {
347 hist << printHeader();
348 }
349 if ( algo_state.iter == 0 ) {
350 hist << " ";
351 hist << std::setw(6) << std::left << algo_state.iter;
352 hist << std::setw(15) << std::left << algo_state.value;
353 hist << std::setw(15) << std::left << algo_state.gnorm;
354 hist << "\n";
355 }
356 else {
357 hist << " ";
358 hist << std::setw(6) << std::left << algo_state.iter;
359 hist << std::setw(15) << std::left << algo_state.value;
360 hist << std::setw(15) << std::left << algo_state.gnorm;
361 hist << std::setw(15) << std::left << algo_state.snorm;
362 hist << std::setw(10) << std::left << algo_state.nfval;
363 hist << std::setw(10) << std::left << algo_state.ngrad;
364 hist << std::setw(10) << std::left << iterKrylov_;
365 hist << std::setw(10) << std::left << flagKrylov_;
366 hist << "\n";
367 }
368 return hist.str();
369 }
370}; // class ProjectedNewtonKrylovStep
371
372} // namespace ROL
373
374#endif
Contains definitions of custom data types in ROL.
Provides the interface to apply upper and lower bound constraints.
void computeProjectedGradient(Vector< Real > &g, const Vector< Real > &x)
Compute projected gradient.
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
Provides definitions for Krylov solvers.
Provides the interface to apply a linear operator.
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 ROL::Ptr< BoundConstraint< Real > > bnd_
void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply linear operator.
HessianPNK(const ROL::Ptr< Objective< Real > > &obj, const ROL::Ptr< BoundConstraint< Real > > &bnd, const ROL::Ptr< Vector< Real > > &x, const ROL::Ptr< Vector< Real > > &g, Real eps=0)
void applyInverse(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply inverse of linear operator.
void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply linear operator.
const ROL::Ptr< BoundConstraint< Real > > bnd_
PrecondPNK(const ROL::Ptr< Secant< Real > > &secant, const ROL::Ptr< BoundConstraint< Real > > &bnd, const ROL::Ptr< Vector< Real > > &x, const ROL::Ptr< Vector< Real > > &g, Real eps=0)
PrecondPNK(const ROL::Ptr< Objective< Real > > &obj, const ROL::Ptr< BoundConstraint< Real > > &bnd, const ROL::Ptr< Vector< Real > > &x, const ROL::Ptr< Vector< Real > > &g, Real eps=0)
Provides the interface to compute optimization steps with projected inexact ProjectedNewton's method ...
std::string print(AlgorithmState< Real > &algo_state, bool print_header=false) const
Print iterate status.
ROL::Ptr< Secant< Real > > secant_
Secant object (used for quasi-Newton)
std::string printName(void) const
Print step name.
bool useSecantPrecond_
Whether or not a secant approximation is used for preconditioning inexact Newton.
void initialize(Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Initialize step with bound constraint.
int flagKrylov_
Termination flag for Krylov method (used for inexact Newton)
bool useProjectedGrad_
Whether or not to use to the projected gradient criticality measure.
void update(Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Update step, if successful.
void compute(Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Compute step.
ProjectedNewtonKrylovStep(ROL::ParameterList &parlist, const ROL::Ptr< Krylov< Real > > &krylov, const ROL::Ptr< Secant< Real > > &secant, const bool computeObj=true)
Constructor.
std::string printHeader(void) const
Print iterate header.
ProjectedNewtonKrylovStep(ROL::ParameterList &parlist, const bool computeObj=true)
Constructor.
int iterKrylov_
Number of Krylov iterations (used for inexact Newton)
ROL::Ptr< Krylov< Real > > krylov_
Krylov solver object (used for inexact Newton)
Provides interface for and implements limited-memory secant operators.
Provides the interface to compute optimization steps.
Definition ROL_Step.hpp:34
virtual void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Initialize step with bound constraint.
Definition ROL_Step.hpp:54
ROL::Ptr< StepState< Real > > getState(void)
Definition ROL_Step.hpp:39
Defines the linear algebra or vector space interface.
virtual Real norm() const =0
Returns where .
virtual void set(const Vector &x)
Set where .
virtual void scale(const Real alpha)=0
Compute where .
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
virtual void plus(const Vector &x)=0
Compute , where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
EKrylov StringToEKrylov(std::string s)
@ DESCENT_NEWTONKRYLOV
ESecant StringToESecant(std::string s)
@ SECANT_USERDEFINED
std::string EDescentToString(EDescent tr)
State for algorithm class. Will be used for restarts.
ROL::Ptr< Vector< Real > > iterateVec