ROL
ROL_ProjectedNewtonKrylovStep.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_PROJECTEDNEWTONKRYLOVSTEP_H
45 #define ROL_PROJECTEDNEWTONKRYLOVSTEP_H
46 
47 #include "ROL_Types.hpp"
48 #include "ROL_Step.hpp"
49 
50 #include "ROL_Secant.hpp"
51 #include "ROL_Krylov.hpp"
52 #include "ROL_LinearOperator.hpp"
53 
54 #include <sstream>
55 #include <iomanip>
56 
64 namespace ROL {
65 
66 template <class Real>
67 class ProjectedNewtonKrylovStep : public Step<Real> {
68 private:
69 
70  ROL::Ptr<Secant<Real> > secant_;
71  ROL::Ptr<Krylov<Real> > krylov_;
72 
75 
76  ROL::Ptr<Vector<Real> > gp_;
77  ROL::Ptr<Vector<Real> > d_;
78 
81  int verbosity_;
82  const bool computeObj_;
83 
86 
87  std::string krylovName_;
88  std::string secantName_;
89 
90  class HessianPNK : public LinearOperator<Real> {
91  private:
92  const ROL::Ptr<Objective<Real> > obj_;
93  const ROL::Ptr<BoundConstraint<Real> > bnd_;
94  const ROL::Ptr<Vector<Real> > x_;
95  const ROL::Ptr<Vector<Real> > g_;
96  ROL::Ptr<Vector<Real> > v_;
97  Real eps_;
98  public:
99  HessianPNK(const ROL::Ptr<Objective<Real> > &obj,
100  const ROL::Ptr<BoundConstraint<Real> > &bnd,
101  const ROL::Ptr<Vector<Real> > &x,
102  const ROL::Ptr<Vector<Real> > &g,
103  Real eps = 0 )
104  : obj_(obj), bnd_(bnd), x_(x), g_(g), eps_(eps) {
105  v_ = x_->clone();
106  }
107  void apply(Vector<Real> &Hv, const Vector<Real> &v, Real &tol) const {
108  v_->set(v);
109  bnd_->pruneActive(*v_,*g_,*x_,eps_);
110  obj_->hessVec(Hv,*v_,*x_,tol);
111  bnd_->pruneActive(Hv,*g_,*x_,eps_);
112  v_->set(v);
113  bnd_->pruneInactive(*v_,*g_,*x_,eps_);
114  Hv.plus(v_->dual());
115  }
116  };
117 
118  class PrecondPNK : public LinearOperator<Real> {
119  private:
120  const ROL::Ptr<Objective<Real> > obj_;
121  const ROL::Ptr<Secant<Real> > secant_;
122  const ROL::Ptr<BoundConstraint<Real> > bnd_;
123  const ROL::Ptr<Vector<Real> > x_;
124  const ROL::Ptr<Vector<Real> > g_;
125  ROL::Ptr<Vector<Real> > v_;
126  Real eps_;
127  const bool useSecant_;
128  public:
129  PrecondPNK(const ROL::Ptr<Objective<Real> > &obj,
130  const ROL::Ptr<BoundConstraint<Real> > &bnd,
131  const ROL::Ptr<Vector<Real> > &x,
132  const ROL::Ptr<Vector<Real> > &g,
133  Real eps = 0 )
134  : obj_(obj), bnd_(bnd), x_(x), g_(g), eps_(eps), useSecant_(false) {
135  v_ = x_->clone();
136  }
137  PrecondPNK(const ROL::Ptr<Secant<Real> > &secant,
138  const ROL::Ptr<BoundConstraint<Real> > &bnd,
139  const ROL::Ptr<Vector<Real> > &x,
140  const ROL::Ptr<Vector<Real> > &g,
141  Real eps = 0 )
142  : secant_(secant), bnd_(bnd), x_(x), g_(g), eps_(eps), useSecant_(true) {
143  v_ = x_->clone();
144  }
145  void apply(Vector<Real> &Hv, const Vector<Real> &v, Real &tol) const {
146  Hv.set(v.dual());
147  }
148  void applyInverse(Vector<Real> &Hv, const Vector<Real> &v, Real &tol) const {
149  v_->set(v);
150  bnd_->pruneActive(*v_,*g_,*x_,eps_);
151  if ( useSecant_ ) {
152  secant_->applyH(Hv,*v_);
153  }
154  else {
155  obj_->precond(Hv,*v_,*x_,tol);
156  }
157  bnd_->pruneActive(Hv,*g_,*x_,eps_);
158  v_->set(v);
159  bnd_->pruneInactive(*v_,*g_,*x_,eps_);
160  Hv.plus(v_->dual());
161  }
162  };
163 
164 public:
165 
167  using Step<Real>::compute;
168  using Step<Real>::update;
169 
177  ProjectedNewtonKrylovStep( ROL::ParameterList &parlist, const bool computeObj = true )
178  : Step<Real>(), secant_(ROL::nullPtr), krylov_(ROL::nullPtr),
179  gp_(ROL::nullPtr), d_(ROL::nullPtr),
181  computeObj_(computeObj), useSecantPrecond_(false) {
182  // Parse ParameterList
183  ROL::ParameterList& Glist = parlist.sublist("General");
184  useSecantPrecond_ = Glist.sublist("Secant").get("Use as Preconditioner", false);
185  useProjectedGrad_ = Glist.get("Projected Gradient Criticality Measure", false);
186  verbosity_ = Glist.get("Print Verbosity",0);
187  // Initialize Krylov object
188  krylovName_ = Glist.sublist("Krylov").get("Type","Conjugate Gradients");
190  krylov_ = KrylovFactory<Real>(parlist);
191  // Initialize secant object
192  secantName_ = Glist.sublist("Secant").get("Type","Limited-Memory BFGS");
193  esec_ = StringToESecant(secantName_);
194  if ( useSecantPrecond_ ) {
195  secant_ = SecantFactory<Real>(parlist);
196  }
197  }
198 
209  ProjectedNewtonKrylovStep(ROL::ParameterList &parlist,
210  const ROL::Ptr<Krylov<Real> > &krylov,
211  const ROL::Ptr<Secant<Real> > &secant,
212  const bool computeObj = true)
213  : Step<Real>(), secant_(secant), krylov_(krylov),
215  gp_(ROL::nullPtr), d_(ROL::nullPtr),
217  computeObj_(computeObj), useSecantPrecond_(false) {
218  // Parse ParameterList
219  ROL::ParameterList& Glist = parlist.sublist("General");
220  useSecantPrecond_ = Glist.sublist("Secant").get("Use as Preconditioner", false);
221  useProjectedGrad_ = Glist.get("Projected Gradient Criticality Measure", false);
222  verbosity_ = Glist.get("Print Verbosity",0);
223  // Initialize secant object
224  if ( useSecantPrecond_ ) {
225  if (secant_ == ROL::nullPtr ) {
226  secantName_ = Glist.sublist("Secant").get("Type","Limited-Memory BFGS");
228  secant_ = SecantFactory<Real>(parlist);
229  }
230  else {
231  secantName_ = Glist.sublist("Secant").get("User Defined Secant Name",
232  "Unspecified User Defined Secant Method");
233  }
234  }
235  // Initialize Krylov object
236  if ( krylov_ == ROL::nullPtr ) {
237  krylovName_ = Glist.sublist("Krylov").get("Type","Conjugate Gradients");
239  krylov_ = KrylovFactory<Real>(parlist);
240  }
241  }
242 
243  void initialize( Vector<Real> &x, const Vector<Real> &s, const Vector<Real> &g,
245  AlgorithmState<Real> &algo_state ) {
246  Step<Real>::initialize(x,s,g,obj,bnd,algo_state);
247  gp_ = g.clone();
248  d_ = s.clone();
249  }
250 
251  void compute( Vector<Real> &s, const Vector<Real> &x,
253  AlgorithmState<Real> &algo_state ) {
254  Real one(1);
255  ROL::Ptr<StepState<Real> > step_state = Step<Real>::getState();
256 
257  // Build Hessian and Preconditioner object
258  ROL::Ptr<Objective<Real> > obj_ptr = ROL::makePtrFromRef(obj);
259  ROL::Ptr<BoundConstraint<Real> > bnd_ptr = ROL::makePtrFromRef(bnd);
260  ROL::Ptr<LinearOperator<Real> > hessian
261  = ROL::makePtr<HessianPNK>(obj_ptr,bnd_ptr,algo_state.iterateVec,
262  step_state->gradientVec,algo_state.gnorm);
263  ROL::Ptr<LinearOperator<Real> > precond;
264  if (useSecantPrecond_) {
265  precond = ROL::makePtr<PrecondPNK>(secant_,bnd_ptr,
266  algo_state.iterateVec,step_state->gradientVec,algo_state.gnorm);
267  }
268  else {
269  precond = ROL::makePtr<PrecondPNK>(obj_ptr,bnd_ptr,
270  algo_state.iterateVec,step_state->gradientVec,algo_state.gnorm);
271  }
272 
273  // Run Krylov method
274  flagKrylov_ = 0;
275  krylov_->run(s,*hessian,*(step_state->gradientVec),*precond,iterKrylov_,flagKrylov_);
276 
277  // Check Krylov flags
278  if ( flagKrylov_ == 2 && iterKrylov_ <= 1 ) {
279  s.set((step_state->gradientVec)->dual());
280  }
281  s.scale(-one);
282  }
283 
284  void update( Vector<Real> &x, const Vector<Real> &s,
286  AlgorithmState<Real> &algo_state ) {
287  Real tol = std::sqrt(ROL_EPSILON<Real>()), one(1);
288  ROL::Ptr<StepState<Real> > step_state = Step<Real>::getState();
289  step_state->SPiter = iterKrylov_;
290  step_state->SPflag = flagKrylov_;
291 
292  // Update iterate and store previous step
293  algo_state.iter++;
294  d_->set(x);
295  x.plus(s);
296  bnd.project(x);
297  (step_state->descentVec)->set(x);
298  (step_state->descentVec)->axpy(-one,*d_);
299  algo_state.snorm = s.norm();
300 
301  // Compute new gradient
302  if ( useSecantPrecond_ ) {
303  gp_->set(*(step_state->gradientVec));
304  }
305  obj.update(x,true,algo_state.iter);
306  if ( computeObj_ ) {
307  algo_state.value = obj.value(x,tol);
308  algo_state.nfval++;
309  }
310  obj.gradient(*(step_state->gradientVec),x,tol);
311  algo_state.ngrad++;
312 
313  // Update Secant Information
314  if ( useSecantPrecond_ ) {
315  secant_->updateStorage(x,*(step_state->gradientVec),*gp_,s,algo_state.snorm,algo_state.iter+1);
316  }
317 
318  // Update algorithm state
319  (algo_state.iterateVec)->set(x);
320  if ( useProjectedGrad_ ) {
321  gp_->set(*(step_state->gradientVec));
322  bnd.computeProjectedGradient( *gp_, x );
323  algo_state.gnorm = gp_->norm();
324  }
325  else {
326  d_->set(x);
327  d_->axpy(-one,(step_state->gradientVec)->dual());
328  bnd.project(*d_);
329  d_->axpy(-one,x);
330  algo_state.gnorm = d_->norm();
331  }
332  }
333 
334  std::string printHeader( void ) const {
335  std::stringstream hist;
336 
337  if( verbosity_>0 ) {
338  hist << std::string(109,'-') << "\n";
340  hist << " status output definitions\n\n";
341  hist << " iter - Number of iterates (steps taken) \n";
342  hist << " value - Objective function value \n";
343  hist << " gnorm - Norm of the gradient\n";
344  hist << " snorm - Norm of the step (update to optimization vector)\n";
345  hist << " #fval - Cumulative number of times the objective function was evaluated\n";
346  hist << " #grad - Number of times the gradient was computed\n";
347  hist << " iterCG - Number of Krylov iterations used to compute search direction\n";
348  hist << " flagCG - Krylov solver flag" << "\n";
349  hist << std::string(109,'-') << "\n";
350  }
351 
352  hist << " ";
353  hist << std::setw(6) << std::left << "iter";
354  hist << std::setw(15) << std::left << "value";
355  hist << std::setw(15) << std::left << "gnorm";
356  hist << std::setw(15) << std::left << "snorm";
357  hist << std::setw(10) << std::left << "#fval";
358  hist << std::setw(10) << std::left << "#grad";
359  hist << std::setw(10) << std::left << "iterCG";
360  hist << std::setw(10) << std::left << "flagCG";
361  hist << "\n";
362  return hist.str();
363  }
364  std::string printName( void ) const {
365  std::stringstream hist;
366  hist << "\n" << EDescentToString(DESCENT_NEWTONKRYLOV);
367  hist << " using " << krylovName_;
368  if ( useSecantPrecond_ ) {
369  hist << " with " << secantName_ << " preconditioning";
370  }
371  hist << "\n";
372  return hist.str();
373  }
374  std::string print( AlgorithmState<Real> &algo_state, bool print_header = false ) const {
375  std::stringstream hist;
376  hist << std::scientific << std::setprecision(6);
377  if ( algo_state.iter == 0 ) {
378  hist << printName();
379  }
380  if ( print_header ) {
381  hist << printHeader();
382  }
383  if ( algo_state.iter == 0 ) {
384  hist << " ";
385  hist << std::setw(6) << std::left << algo_state.iter;
386  hist << std::setw(15) << std::left << algo_state.value;
387  hist << std::setw(15) << std::left << algo_state.gnorm;
388  hist << "\n";
389  }
390  else {
391  hist << " ";
392  hist << std::setw(6) << std::left << algo_state.iter;
393  hist << std::setw(15) << std::left << algo_state.value;
394  hist << std::setw(15) << std::left << algo_state.gnorm;
395  hist << std::setw(15) << std::left << algo_state.snorm;
396  hist << std::setw(10) << std::left << algo_state.nfval;
397  hist << std::setw(10) << std::left << algo_state.ngrad;
398  hist << std::setw(10) << std::left << iterKrylov_;
399  hist << std::setw(10) << std::left << flagKrylov_;
400  hist << "\n";
401  }
402  return hist.str();
403  }
404 }; // class ProjectedNewtonKrylovStep
405 
406 } // namespace ROL
407 
408 #endif
Provides the interface to evaluate objective functions.
void applyInverse(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply inverse of linear operator.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: ROL_Vector.hpp:226
bool useSecantPrecond_
Whether or not a secant approximation is used for preconditioning inexact Newton. ...
ROL::Ptr< Secant< Real > > secant_
Secant object (used for quasi-Newton)
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
void update(Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Update step, if successful.
virtual void plus(const Vector &x)=0
Compute , where .
const ROL::Ptr< BoundConstraint< Real > > bnd_
void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply linear operator.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
Provides the interface to compute optimization steps.
Definition: ROL_Step.hpp:69
Contains definitions of custom data types in ROL.
void compute(Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Compute step.
void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply linear operator.
ESecant StringToESecant(std::string s)
Definition: ROL_Types.hpp:541
std::string EDescentToString(EDescent tr)
Definition: ROL_Types.hpp:418
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
Provides the interface to compute optimization steps with projected inexact ProjectedNewton&#39;s method ...
ProjectedNewtonKrylovStep(ROL::ParameterList &parlist, const bool computeObj=true)
Constructor.
EKrylov
Enumeration of Krylov methods.
Definition: ROL_Types.hpp:557
std::string printName(void) const
Print step name.
EKrylov StringToEKrylov(std::string s)
Definition: ROL_Types.hpp:612
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)
ProjectedNewtonKrylovStep(ROL::ParameterList &parlist, const ROL::Ptr< Krylov< Real > > &krylov, const ROL::Ptr< Secant< Real > > &secant, const bool computeObj=true)
Constructor.
State for algorithm class. Will be used for restarts.
Definition: ROL_Types.hpp:143
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)
std::string printHeader(void) const
Print iterate header.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
bool useProjectedGrad_
Whether or not to use to the projected gradient criticality measure.
ESecant
Enumeration of secant update algorithms.
Definition: ROL_Types.hpp:484
ROL::Ptr< StepState< Real > > getState(void)
Definition: ROL_Step.hpp:74
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)
Provides interface for and implements limited-memory secant operators.
Definition: ROL_Secant.hpp:70
ROL::Ptr< Vector< Real > > iterateVec
Definition: ROL_Types.hpp:157
Provides definitions for Krylov solvers.
Definition: ROL_Krylov.hpp:58
Provides the interface to apply a linear operator.
Provides the interface to apply upper and lower bound constraints.
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 computeProjectedGradient(Vector< Real > &g, const Vector< Real > &x)
Compute projected gradient.
int iterKrylov_
Number of Krylov iterations (used for inexact Newton)
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:89
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
const ROL::Ptr< BoundConstraint< Real > > bnd_
virtual Real norm() const =0
Returns where .
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
std::string print(AlgorithmState< Real > &algo_state, bool print_header=false) const
Print iterate status.
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
ROL::Ptr< Krylov< Real > > krylov_
Krylov solver object (used for inexact Newton)