44 #ifndef ROL_PROJECTEDNEWTONKRYLOVSTEP_H
45 #define ROL_PROJECTEDNEWTONKRYLOVSTEP_H
76 ROL::Ptr<Vector<Real> >
gp_;
77 ROL::Ptr<Vector<Real> >
d_;
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_;
120 const ROL::Ptr<Objective<Real> >
obj_;
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_;
179 gp_(ROL::nullPtr),
d_(ROL::nullPtr),
183 ROL::ParameterList& Glist = parlist.sublist(
"General");
188 krylovName_ = Glist.sublist(
"Krylov").get(
"Type",
"Conjugate Gradients");
190 krylov_ = KrylovFactory<Real>(parlist);
192 secantName_ = Glist.sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS");
195 secant_ = SecantFactory<Real>(parlist);
212 const bool computeObj =
true)
215 gp_(ROL::nullPtr),
d_(ROL::nullPtr),
219 ROL::ParameterList& Glist = parlist.sublist(
"General");
225 if (
secant_ == ROL::nullPtr ) {
226 secantName_ = Glist.sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS");
228 secant_ = SecantFactory<Real>(parlist);
231 secantName_ = Glist.sublist(
"Secant").get(
"User Defined Secant Name",
232 "Unspecified User Defined Secant Method");
236 if (
krylov_ == ROL::nullPtr ) {
237 krylovName_ = Glist.sublist(
"Krylov").get(
"Type",
"Conjugate Gradients");
239 krylov_ = KrylovFactory<Real>(parlist);
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;
265 precond = ROL::makePtr<PrecondPNK>(
secant_,bnd_ptr,
269 precond = ROL::makePtr<PrecondPNK>(obj_ptr,bnd_ptr,
279 s.
set((step_state->gradientVec)->dual());
287 Real tol = std::sqrt(ROL_EPSILON<Real>()), one(1);
297 (step_state->descentVec)->set(x);
298 (step_state->descentVec)->axpy(-one,*
d_);
303 gp_->set(*(step_state->gradientVec));
310 obj.
gradient(*(step_state->gradientVec),x,tol);
315 secant_->updateStorage(x,*(step_state->gradientVec),*
gp_,s,algo_state.
snorm,algo_state.
iter+1);
321 gp_->set(*(step_state->gradientVec));
327 d_->axpy(-one,(step_state->gradientVec)->dual());
335 std::stringstream hist;
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";
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";
365 std::stringstream hist;
369 hist <<
" with " <<
secantName_ <<
" preconditioning";
375 std::stringstream hist;
376 hist << std::scientific << std::setprecision(6);
377 if ( algo_state.
iter == 0 ) {
380 if ( print_header ) {
383 if ( algo_state.
iter == 0 ) {
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;
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;
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...
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.
const ROL::Ptr< Vector< Real > > g_
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
Provides the interface to compute optimization steps.
const ROL::Ptr< Objective< Real > > obj_
ROL::Ptr< Vector< Real > > v_
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)
std::string EDescentToString(EDescent tr)
Defines the linear algebra or vector space interface.
Provides the interface to compute optimization steps with projected inexact ProjectedNewton's method ...
const ROL::Ptr< Vector< Real > > g_
ProjectedNewtonKrylovStep(ROL::ParameterList &parlist, const bool computeObj=true)
Constructor.
EKrylov
Enumeration of Krylov methods.
std::string printName(void) const
Print step name.
EKrylov StringToEKrylov(std::string s)
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.
ROL::Ptr< Vector< Real > > d_
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.
const ROL::Ptr< Vector< Real > > x_
ROL::Ptr< StepState< Real > > getState(void)
const ROL::Ptr< Secant< Real > > secant_
ROL::Ptr< Vector< Real > > v_
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.
ROL::Ptr< Vector< Real > > iterateVec
Provides definitions for Krylov solvers.
Provides the interface to apply a linear operator.
ROL::Ptr< Vector< Real > > gp_
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.
const ROL::Ptr< Vector< Real > > x_
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.
virtual void set(const Vector &x)
Set where .
const ROL::Ptr< BoundConstraint< Real > > bnd_
virtual Real norm() const =0
Returns where .
const ROL::Ptr< Objective< Real > > obj_
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)
int verbosity_
Verbosity level.