10#ifndef ROL_TYPEB_LSecantBALGORITHM_DEF_HPP
11#define ROL_TYPEB_LSecantBALGORITHM_DEF_HPP
16template<
typename Real>
23 ParameterList &trlist = list.sublist(
"Step").sublist(
"Line Search");
24 state_->searchSize =
static_cast<Real
>(1);
26 maxit_ = list.sublist(
"General").sublist(
"Krylov").get(
"Iteration Limit", 20);
27 tol1_ = list.sublist(
"General").sublist(
"Krylov").get(
"Absolute Tolerance", 1e-4);
28 tol2_ = list.sublist(
"General").sublist(
"Krylov").get(
"Relative Tolerance", 1e-2);
30 ROL::ParameterList &lmlist = trlist.sublist(
"Quasi-Newton").sublist(
"L-Secant-B");
31 mu0_ = lmlist.get(
"Sufficient Decrease Parameter", 1e-2);
32 spexp_ = lmlist.get(
"Relative Tolerance Exponent", 1.0);
33 spexp_ = std::max(
static_cast<Real
>(1),std::min(
spexp_,
static_cast<Real
>(2)));
34 redlim_ = lmlist.sublist(
"Cauchy Point").get(
"Maximum Number of Reduction Steps", 10);
35 explim_ = lmlist.sublist(
"Cauchy Point").get(
"Maximum Number of Expansion Steps", 10);
36 alpha_ = lmlist.sublist(
"Cauchy Point").get(
"Initial Step Size", 1.0);
37 normAlpha_ = lmlist.sublist(
"Cauchy Point").get(
"Normalize Initial Step Size",
false);
38 interpf_ = lmlist.sublist(
"Cauchy Point").get(
"Reduction Rate", 0.1);
39 extrapf_ = lmlist.sublist(
"Cauchy Point").get(
"Expansion Rate", 10.0);
40 qtol_ = lmlist.sublist(
"Cauchy Point").get(
"Decrease Tolerance", 1e-8);
41 interpfPS_ = trlist.sublist(
"Line-Search Method").get(
"Backtracking Rate", 0.5);
43 verbosity_ = list.sublist(
"General").get(
"Output Level",0);
49 if (secant == nullPtr) {
50 std::string secantType = list.sublist(
"General").sublist(
"Secant").get(
"Type",
"Limited-Memory Secant");
56template<
typename Real>
61 std::ostream &outStream) {
64 if (
proj_ == nullPtr) {
65 proj_ = makePtr<PolyhedralProjection<Real>>(makePtrFromRef(bnd));
73 state_->iterateVec->set(x);
80 state_->stepVec->axpy(-one,x);
89 rcon_ = makePtr<ReducedLinearConstraint<Real>>(
proj_->getLinearConstraint(),
92 ns_ = makePtr<NullSpaceOperator<Real>>(
rcon_,x,
93 *
proj_->getResidual());
97template<
typename Real>
102 std::ostream &outStream ) {
103 const Real
zero(0), one(1);
105 Real gfnorm(0), gs(0), ftrial(0), q(0), tol(0), stol(0), snorm(0);
108 Ptr<Vector<Real>> s = x.
clone();
109 Ptr<Vector<Real>> gmod = g.
clone(), gfree = g.
clone();
110 Ptr<Vector<Real>> pwa1 = x.
clone(), pwa2 = x.
clone(), pwa3 = x.
clone();
111 Ptr<Vector<Real>> dwa1 = g.
clone(), dwa2 = g.
clone(), dwa3 = g.
clone();
123 state_->gradientVec->dual(),
132 gmod->plus(*
state_->gradientVec);
135 pwa1->set(gfree->dual());
137 gfree->set(pwa1->dual());
139 gfnorm = gfree->norm();
142 gfnorm = pwa1->norm();
146 gfnorm = gfree->norm();
150 outStream <<
" Norm of free gradient components: " << gfnorm << std::endl;
151 outStream << std::endl;
161 *pwa1, *dwa1, *pwa2, *dwa2, *pwa3, *dwa3,
164 outStream <<
" Computation of CG step" << std::endl;
165 outStream <<
" CG step length: " << snorm << std::endl;
166 outStream <<
" Number of CG iterations: " <<
SPiter_ << std::endl;
167 outStream <<
" CG flag: " <<
SPflag_ << std::endl;
168 outStream << std::endl;
181 state_->value,gs,obj,*pwa1,outStream);
188 dwa1->set(*
state_->gradientVec);
202template<
typename Real>
205 std::ostream &outStream)
const {
213template<
typename Real>
221 std::ostream &outStream) {
222 const Real half(0.5);
224 Real gs(0), snorm(0);
226 snorm =
dgpstep(s,g,x,-alpha,outStream);
229 q = half * s.
apply(dwa) + gs;
230 interp = (q >
mu0_*gs);
237 snorm =
dgpstep(s,g,x,-alpha,outStream);
240 q = half * s.
apply(dwa) + gs;
252 snorm =
dgpstep(s,g,x,-alpha,outStream);
256 q = half * s.
apply(dwa) + gs;
257 if (q <=
mu0_*gs && std::abs(q-qs) >
qtol_*std::abs(qs)) {
275 snorm =
dgpstep(s,g,x,-alpha,outStream);
278 outStream <<
" Cauchy point" << std::endl;
279 outStream <<
" Step length (alpha): " << alpha << std::endl;
280 outStream <<
" Step length (alpha*g): " << snorm << std::endl;
281 outStream <<
" Model decrease: " << -q << std::endl;
283 outStream <<
" Number of extrapolation steps: " << cnt << std::endl;
289template<
typename Real>
305 if (fnew <= fold +
mu0_*beta*gs) search =
false;
312 outStream << std::endl;
313 outStream <<
" Line search" << std::endl;
314 outStream <<
" Step length (beta): " << beta << std::endl;
315 outStream <<
" Step length (beta*s): " << snorm << std::endl;
316 outStream <<
" New objective value: " << fnew << std::endl;
317 outStream <<
" Old objective value: " << fold << std::endl;
318 outStream <<
" Descent verification (gs): " << gs << std::endl;
319 outStream <<
" Number of steps: " << nsteps << std::endl;
324template<
typename Real>
329 const Real tol,
const Real stol,
const int itermax,
332 std::ostream &outStream)
const {
339 const Real
zero(0), half(0.5), one(1), two(2);
340 Real rho(0), kappa(0), beta(0), sigma(0), alpha(0), alphatmp(0);
341 Real rtr(0), tnorm(0), sMs(0), pMp(0), sMp(0);
355 for (iter = 0; iter < itermax; ++iter) {
365 tolB = minIntSize * alpha;
366 while (alpha-sigma > tolB) {
367 alphatmp = sigma + half * (alpha - sigma);
370 else sigma = alphatmp;
374 sMs = sMs + two*sigma*sMp + sigma*sigma*pMp;
386 if (rtr <= stol*stol || tnorm <= tol) {
387 sMs = sMs + two*alpha*sMp + alpha*alpha*pMp;
399 sMs = sMs + two*alpha*sMp + alpha*alpha*pMp;
400 sMp = beta*(sMp + alpha*pMp);
401 pMp = (!
hasEcon_ ? rho : p.
dot(p)) + beta*beta*pMp;
404 if (iter == itermax) iflag = 1;
405 if (iflag != 1) iter++;
406 return std::sqrt(sMs);
409template<
typename Real>
426template<
typename Real>
445 rcon_->setX(makePtrFromRef(x));
448 ns_->apply(hv,pwa,tol);
452template<
typename Real>
454 std::ios_base::fmtflags osFlags(os.flags());
456 os << std::string(114,
'-') << std::endl;
457 os <<
" L-Secant-B line search method status output definitions" << std::endl << std::endl;
458 os <<
" iter - Number of iterates (steps taken)" << std::endl;
459 os <<
" value - Objective function value" << std::endl;
460 os <<
" gnorm - Norm of the gradient" << std::endl;
461 os <<
" snorm - Norm of the step (update to optimization vector)" << std::endl;
462 os <<
" LSpar - Line-Search parameter" << std::endl;
463 os <<
" #fval - Number of times the objective function was evaluated" << std::endl;
464 os <<
" #grad - Number of times the gradient was computed" << std::endl;
465 os <<
" #proj - Number of times the projection was applied" << std::endl;
466 os <<
" iterCG - Number of Truncated CG iterations" << std::endl << std::endl;
467 os <<
" flagGC - Trust-Region Truncated CG flag" << std::endl;
468 os <<
" 0 - Converged" << std::endl;
469 os <<
" 1 - Iteration Limit Exceeded" << std::endl;
470 os <<
" 2 - Bounds Exceeded" << std::endl;
471 os << std::string(114,
'-') << std::endl;
474 os << std::setw(6) << std::left <<
"iter";
475 os << std::setw(15) << std::left <<
"value";
476 os << std::setw(15) << std::left <<
"gnorm";
477 os << std::setw(15) << std::left <<
"snorm";
478 os << std::setw(15) << std::left <<
"LSpar";
479 os << std::setw(10) << std::left <<
"#fval";
480 os << std::setw(10) << std::left <<
"#grad";
481 os << std::setw(10) << std::left <<
"#proj";
482 os << std::setw(10) << std::left <<
"iterCG";
483 os << std::setw(10) << std::left <<
"flagCG";
488template<
typename Real>
490 std::ios_base::fmtflags osFlags(os.flags());
491 os << std::endl <<
"L-Secant-B Line-Search Method (Type B, Bound Constraints)" << std::endl;
495template<
typename Real>
497 std::ios_base::fmtflags osFlags(os.flags());
498 os << std::scientific << std::setprecision(6);
501 if (
state_->iter == 0 ) {
503 os << std::setw(6) << std::left <<
state_->iter;
504 os << std::setw(15) << std::left <<
state_->value;
505 os << std::setw(15) << std::left <<
state_->gnorm;
506 os << std::setw(15) << std::left <<
"---";
507 os << std::setw(15) << std::left <<
state_->searchSize;
508 os << std::setw(10) << std::left <<
state_->nfval;
509 os << std::setw(10) << std::left <<
state_->ngrad;
510 os << std::setw(10) << std::left <<
state_->nproj;
511 os << std::setw(10) << std::left <<
"---";
512 os << std::setw(10) << std::left <<
"---";
517 os << std::setw(6) << std::left <<
state_->iter;
518 os << std::setw(15) << std::left <<
state_->value;
519 os << std::setw(15) << std::left <<
state_->gnorm;
520 os << std::setw(15) << std::left <<
state_->snorm;
521 os << std::setw(15) << std::left <<
state_->searchSize;
522 os << std::setw(10) << std::left <<
state_->nfval;
523 os << std::setw(10) << std::left <<
state_->ngrad;
524 os << std::setw(10) << std::left <<
state_->nproj;
525 os << std::setw(10) << std::left <<
SPiter_;
526 os << std::setw(10) << std::left <<
SPflag_;
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0 zero)()
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.
Provides the interface to apply upper and lower bound constraints.
virtual bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the -active set.
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.
virtual void applyB(Vector< Real > &Bv, const Vector< Real > &v) const =0
virtual void applyH(Vector< Real > &Hv, const Vector< Real > &v) const =0
Provides an interface to check status of optimization algorithms.
Ptr< PolyhedralProjection< Real > > proj_
void initialize(const Vector< Real > &x, const Vector< Real > &g)
Real optimalityCriterion(const Vector< Real > &x, const Vector< Real > &g, Vector< Real > &primal, std::ostream &outStream=std::cout) const
virtual void writeExitStatus(std::ostream &os) const
const Ptr< AlgorithmState< Real > > state_
const Ptr< CombinedStatusTest< Real > > status_
Real alpha_
Initial Cauchy point step length (default: 1.0).
Ptr< ReducedLinearConstraint< Real > > rcon_
Equality constraint restricted to current active variables.
int SPiter_
Subproblem solver iteration count.
ESecant esec_
Secant type (default: Limited-Memory BFGS).
bool normAlpha_
Normalize initial Cauchy point step length (default: false).
Ptr< NullSpaceOperator< Real > > ns_
Null space projection onto reduced equality constraint Jacobian.
bool writeHeader_
Flag to write header at every iteration.
void writeOutput(std::ostream &os, const bool write_header=false) const override
Print iterate status.
void applyFreeHessian(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Secant< Real > &secant, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &pwa) const
int redlim_
Maximum number of Cauchy point reduction steps (default: 10).
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...
Real interpfPS_
Backtracking rate for projected search (default: 0.5).
Real interpf_
Backtracking rate for Cauchy point computation (default: 1e-1).
Real spexp_
Relative tolerance exponent for subproblem solve (default: 1, range: [1,2]).
Real dcauchy(Vector< Real > &s, Real &alpha, Real &q, const Vector< Real > &x, const Vector< Real > &g, Secant< Real > &secant, Vector< Real > &dwa, Vector< Real > &dwa1, std::ostream &outStream=std::cout)
Real dpcg(Vector< Real > &w, int &iflag, int &iter, const Vector< Real > &g, const Vector< Real > &x, Secant< Real > &secant, BoundConstraint< Real > &bnd, const Real tol, const Real stol, const int itermax, Vector< Real > &p, Vector< Real > &q, Vector< Real > &r, Vector< Real > &t, Vector< Real > &pwa, Vector< Real > &dwa, std::ostream &outStream=std::cout) const
void applyFreePrecond(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Secant< Real > &secant, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &dwa, Vector< Real > &pwa) const
Ptr< Secant< Real > > secant_
Secant object.
Real dgpstep(Vector< Real > &s, const Vector< Real > &w, const Vector< Real > &x, const Real alpha, std::ostream &outStream=std::cout) const
Real qtol_
Relative tolerance for computed decrease in Cauchy point computation (default: 1-8).
int maxit_
Maximum number of CG iterations (default: 20).
void writeHeader(std::ostream &os) const override
Print iterate header.
Real extrapf_
Extrapolation rate for Cauchy point computation (default: 1e1).
int SPflag_
Subproblem solver termination flag.
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout)
LSecantBAlgorithm(ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr)
Real dsrch(Vector< Real > &x, Vector< Real > &s, Real &fnew, Real &beta, Real fold, Real gs, Objective< Real > &obj, Vector< Real > &pwa, std::ostream &outStream=std::cout)
Real tol2_
Relative tolerance for truncated CG (default: 1e-2).
Real tol1_
Absolute tolerance for truncated CG (default: 1e-4).
bool useSecantPrecond_
Flag to use secant as a preconditioner (default: false).
bool hasEcon_
Flag signifies if equality constraints exist.
Real mu0_
Sufficient decrease parameter (default: 1e-2).
void writeName(std::ostream &os) const override
Print step name.
unsigned verbosity_
Output level (default: 0).
bool useSecantHessVec_
Flag to use secant as Hessian (default: false).
int explim_
Maximum number of Cauchy point expansion steps (default: 10).
Defines the linear algebra or vector space interface.
virtual Real apply(const Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
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 void zero()
Set to zero vector.
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 .
virtual Real dot(const Vector &x) const =0
Compute where .
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
ESecant StringToESecant(std::string s)
Real ROL_OVERFLOW(void)
Platform-dependent maximum double.
ROL::Ptr< Secant< Real > > SecantFactory(ROL::ParameterList &parlist, ESecantMode mode=SECANTMODE_BOTH)