10#ifndef ROL_TYPEB_QUASINEWTONALGORITHM_DEF_HPP
11#define ROL_TYPEB_QUASINEWTONALGORITHM_DEF_HPP
24template<
typename Real>
33 ParameterList &lslist = list.sublist(
"Step").sublist(
"Line Search");
34 maxit_ = lslist.get(
"Function Evaluation Limit", 20);
35 c1_ = lslist.get(
"Sufficient Decrease Tolerance", 1e-4);
36 rhodec_ = lslist.sublist(
"Line-Search Method").get(
"Backtracking Rate", 0.5);
37 sigma1_ = lslist.sublist(
"PQN").get(
"Lower Step Size Safeguard", 0.1);
38 sigma2_ = lslist.sublist(
"PQN").get(
"Upper Step Size Safeguard", 0.9);
39 algoName_ = lslist.sublist(
"PQN").get(
"Subproblem Solver",
"Spectral Gradient");
40 int sp_maxit = lslist.sublist(
"PQN").get(
"Subproblem Iteration Limit", 1000);
41 sp_tol1_ = lslist.sublist(
"PQN").get(
"Subproblem Absolute Tolerance", 1e-4);
42 sp_tol2_ = lslist.sublist(
"PQN").get(
"Subproblem Relative Tolerance", 1e-2);
43 Real opt_tol = lslist.sublist(
"Status Test").get(
"Gradient Tolerance", 1e-8);
45 verbosity_ = list.sublist(
"General").get(
"Output Level", 0);
48 list_.sublist(
"Status Test").set(
"Iteration Limit", sp_maxit);
52 secantName_ = list.sublist(
"General").sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS");
57 secantName_ = list.sublist(
"General").sublist(
"Secant").get(
"User Defined Secant Name",
58 "Unspecified User Defined Secant Method");
63template<
typename Real>
68 std::ostream &outStream) {
70 if (
proj_ == nullPtr) {
71 proj_ = makePtr<PolyhedralProjection<Real>>(makePtrFromRef(bnd));
79 state_->iterateVec->set(x);
86 state_->stepVec->axpy(-one,x);
91template<
typename Real>
96 std::ostream &outStream ) {
97 const Real half(0.5), one(1);
101 Real ftrial(0), gs(0), alphaTmp(0), tol(std::sqrt(
ROL_EPSILON<Real>())), gtol(1);
103 Ptr<TypeB::Algorithm<Real>> algo;
104 Ptr<PQNObjective<Real>> qobj = makePtr<PQNObjective<Real>>(
secant_,x,g);
105 Ptr<Problem<Real>> problem = makePtr<Problem<Real>>(qobj,xs);
106 problem->addBoundConstraint(makePtrFromRef(bnd));
108 problem->addLinearConstraint(
"LEC",
proj_->getLinearConstraint(),
109 proj_->getMultiplier(),
110 proj_->getResidual());
111 problem->setProjectionAlgorithm(
list_);
113 problem->finalize(
false,
verbosity_>2,outStream);
119 gp->set(
state_->gradientVec->dual());
122 qobj->setAnchor(x,*
state_->gradientVec);
123 xs->set(x); xs->axpy(-one,*gp);
proj_->project(*xs,outStream);
state_->nproj++;
125 list_.sublist(
"Status Test").set(
"Gradient Tolerance",gtol);
126 if (
algoName_ ==
"Trust Region") algo = makePtr<TypeB::LinMoreAlgorithm<Real>>(
list_);
127 else if (
algoName_ ==
"Line Search") algo = makePtr<TypeB::GradientAlgorithm<Real>>(
list_);
128 else if (
algoName_ ==
"Primal Dual Active Set") algo = makePtr<TypeB::PrimalDualActiveSetAlgorithm<Real>>(
list_);
129 else if (
algoName_ ==
"Moreau-Yosida") algo = makePtr<TypeB::MoreauYosidaAlgorithm<Real>>(
list_);
130 else if (
algoName_ ==
"Interior Point") algo = makePtr<TypeB::InteriorPointAlgorithm<Real>>(
list_);
131 else algo = makePtr<TypeB::SpectralGradientAlgorithm<Real>>(
list_);
132 algo->run(*problem,outStream);
133 s->set(*xs); s->axpy(-one,x);
135 state_->nproj += staticPtrCast<const TypeB::AlgorithmState<Real>>(algo->getState())->nproj;
143 gs =
state_->gradientVec->apply(*s);
145 outStream <<
" In TypeB::QuasiNewtonAlgorithm: Line Search" << std::endl;
146 outStream <<
" Step size: " <<
state_->searchSize << std::endl;
147 outStream <<
" Trial objective value: " << ftrial << std::endl;
148 outStream <<
" Computed reduction: " <<
state_->value-ftrial << std::endl;
149 outStream <<
" Dot product of gradient and step: " << gs << std::endl;
150 outStream <<
" Sufficient decrease bound: " << -gs*
state_->searchSize*
c1_ << std::endl;
151 outStream <<
" Number of function evaluations: " <<
ls_nfval_ << std::endl;
154 alphaTmp = -half*
state_->searchSize*
state_->searchSize*gs
164 outStream << std::endl;
165 outStream <<
" Step size: " <<
state_->searchSize << std::endl;
166 outStream <<
" Trial objective value: " << ftrial << std::endl;
167 outStream <<
" Computed reduction: " <<
state_->value-ftrial << std::endl;
168 outStream <<
" Dot product of gradient and step: " << gs << std::endl;
169 outStream <<
" Sufficient decrease bound: " << -gs*
state_->searchSize*
c1_ << std::endl;
170 outStream <<
" Number of function evaluations: " <<
ls_nfval_ << std::endl;
181 state_->iterateVec->set(x);
187 gold->set(*
state_->gradientVec);
189 gp->set(
state_->gradientVec->dual());
192 s->set(x); s->axpy(-one,*gp);
195 state_->gnorm = s->norm();
206template<
typename Real>
208 std::ios_base::fmtflags osFlags(os.flags());
210 os << std::string(114,
'-') << std::endl;
211 os <<
"Line-Search Projected Quasi-Newton with " <<
secantName_ <<
" Hessian approximation";
212 os <<
" status output definitions" << std::endl << std::endl;
213 os <<
" iter - Number of iterates (steps taken)" << std::endl;
214 os <<
" value - Objective function value" << std::endl;
215 os <<
" gnorm - Norm of the gradient" << std::endl;
216 os <<
" snorm - Norm of the step (update to optimization vector)" << std::endl;
217 os <<
" alpha - Line search step length" << std::endl;
218 os <<
" #fval - Cumulative number of times the objective function was evaluated" << std::endl;
219 os <<
" #grad - Cumulative number of times the gradient was computed" << std::endl;
220 os <<
" #proj - Cumulative number of times the projection was computed" << std::endl;
221 os <<
" ls_#fval - Number of the times the objective function was evaluated during the line search" << std::endl;
222 os <<
" sp_iter - Number iterations to compute quasi-Newton step" << std::endl;
223 os << std::string(114,
'-') << std::endl;
227 os << std::setw(6) << std::left <<
"iter";
228 os << std::setw(15) << std::left <<
"value";
229 os << std::setw(15) << std::left <<
"gnorm";
230 os << std::setw(15) << std::left <<
"snorm";
231 os << std::setw(15) << std::left <<
"alpha";
232 os << std::setw(10) << std::left <<
"#fval";
233 os << std::setw(10) << std::left <<
"#grad";
234 os << std::setw(10) << std::left <<
"#proj";
235 os << std::setw(10) << std::left <<
"#ls_fval";
236 os << std::setw(10) << std::left <<
"sp_iter";
241template<
typename Real>
243 std::ios_base::fmtflags osFlags(os.flags());
244 os << std::endl <<
"Line-Search Projected Quasi-Newton (Type B, Bound Constraints)" << std::endl;
248template<
typename Real>
250 std::ios_base::fmtflags osFlags(os.flags());
251 os << std::scientific << std::setprecision(6);
254 if (
state_->iter == 0 ) {
256 os << std::setw(6) << std::left <<
state_->iter;
257 os << std::setw(15) << std::left <<
state_->value;
258 os << std::setw(15) << std::left <<
state_->gnorm;
259 os << std::setw(15) << std::left <<
"---";
260 os << std::setw(15) << std::left <<
"---";
261 os << std::setw(10) << std::left <<
state_->nfval;
262 os << std::setw(10) << std::left <<
state_->ngrad;
263 os << std::setw(10) << std::left <<
state_->nproj;
264 os << std::setw(10) << std::left <<
"---";
265 os << std::setw(10) << std::left <<
"---";
270 os << std::setw(6) << std::left <<
state_->iter;
271 os << std::setw(15) << std::left <<
state_->value;
272 os << std::setw(15) << std::left <<
state_->gnorm;
273 os << std::setw(15) << std::left <<
state_->snorm;
274 os << std::setw(15) << std::left <<
state_->searchSize;
275 os << std::setw(10) << std::left <<
state_->nfval;
276 os << std::setw(10) << std::left <<
state_->ngrad;
277 os << std::setw(10) << std::left <<
state_->nproj;
278 os << std::setw(10) << std::left <<
ls_nfval_;
279 os << std::setw(10) << std::left <<
spgIter_;
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.
Provides the interface to apply upper and lower bound constraints.
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.
Provides interface for and implements limited-memory secant operators.
Provides an interface to check status of optimization algorithms.
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_
Real sigma2_
Upper safeguard for quadratic line search (default: 0.9).
Real c1_
Sufficient Decrease Parameter (default: 1e-4).
int maxit_
Maximum number of line search steps (default: 20).
Real sigma1_
Lower safeguard for quadratic line search (default: 0.1).
Ptr< Secant< Real > > secant_
Secant object (used for quasi-Newton).
void writeOutput(std::ostream &os, const bool write_header=false) const override
Print iterate status.
void writeHeader(std::ostream &os) const override
Print iterate header.
void writeName(std::ostream &os) const override
Print step name.
Real rhodec_
Backtracking rate (default: 0.5).
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...
QuasiNewtonAlgorithm(ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr)
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout)
std::string secantName_
Secant name.
ESecant esec_
Secant type.
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 .
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
ESecant StringToESecant(std::string s)
ROL::Ptr< Secant< Real > > SecantFactory(ROL::ParameterList &parlist, ESecantMode mode=SECANTMODE_BOTH)