44#ifndef ROL_TYPEB_LINMOREALGORITHM_DEF_HPP
45#define ROL_TYPEB_LINMOREALGORITHM_DEF_HPP
50template<
typename Real>
57 ParameterList &trlist = list.sublist(
"Step").sublist(
"Trust Region");
59 state_->searchSize = trlist.get(
"Initial Radius", -1.0);
61 eta0_ = trlist.get(
"Step Acceptance Threshold", 0.05);
62 eta1_ = trlist.get(
"Radius Shrinking Threshold", 0.05);
63 eta2_ = trlist.get(
"Radius Growing Threshold", 0.9);
64 gamma0_ = trlist.get(
"Radius Shrinking Rate (Negative rho)", 0.0625);
65 gamma1_ = trlist.get(
"Radius Shrinking Rate (Positive rho)", 0.25);
66 gamma2_ = trlist.get(
"Radius Growing Rate", 2.5);
67 TRsafe_ = trlist.get(
"Safeguard Size", 100.0);
69 interpRad_ = trlist.get(
"Use Radius Interpolation",
false);
71 storageNM_ = trlist.get(
"Nonmonotone Storage Size", 0);
74 maxit_ = list.sublist(
"General").sublist(
"Krylov").get(
"Iteration Limit", 20);
75 tol1_ = list.sublist(
"General").sublist(
"Krylov").get(
"Absolute Tolerance", 1e-4);
76 tol2_ = list.sublist(
"General").sublist(
"Krylov").get(
"Relative Tolerance", 1e-2);
78 ROL::ParameterList &lmlist = trlist.sublist(
"Lin-More");
79 minit_ = lmlist.get(
"Maximum Number of Minor Iterations", 10);
80 mu0_ = lmlist.get(
"Sufficient Decrease Parameter", 1e-2);
81 spexp_ = lmlist.get(
"Relative Tolerance Exponent", 1.0);
82 spexp_ = std::max(
static_cast<Real
>(1),std::min(
spexp_,
static_cast<Real
>(2)));
83 redlim_ = lmlist.sublist(
"Cauchy Point").get(
"Maximum Number of Reduction Steps", 10);
84 explim_ = lmlist.sublist(
"Cauchy Point").get(
"Maximum Number of Expansion Steps", 10);
85 alpha_ = lmlist.sublist(
"Cauchy Point").get(
"Initial Step Size", 1.0);
86 normAlpha_ = lmlist.sublist(
"Cauchy Point").get(
"Normalize Initial Step Size",
false);
87 interpf_ = lmlist.sublist(
"Cauchy Point").get(
"Reduction Rate", 0.1);
88 extrapf_ = lmlist.sublist(
"Cauchy Point").get(
"Expansion Rate", 10.0);
89 qtol_ = lmlist.sublist(
"Cauchy Point").get(
"Decrease Tolerance", 1e-8);
90 interpfPS_ = lmlist.sublist(
"Projected Search").get(
"Backtracking Rate", 0.5);
91 pslim_ = lmlist.sublist(
"Projected Search").get(
"Maximum Number of Steps", 20);
93 verbosity_ = list.sublist(
"General").get(
"Output Level",0);
96 useSecantPrecond_ = list.sublist(
"General").sublist(
"Secant").get(
"Use as Preconditioner",
false);
97 useSecantHessVec_ = list.sublist(
"General").sublist(
"Secant").get(
"Use as Hessian",
false);
102 model_ = makePtr<TrustRegionModel_U<Real>>(list,secant,mode);
103 if (secant == nullPtr) {
104 esec_ =
StringToESecant(list.sublist(
"General").sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS"));
108template<
typename Real>
113 std::ostream &outStream) {
116 if (
proj_ == nullPtr) {
117 proj_ = makePtr<PolyhedralProjection<Real>>(makePtrFromRef(bnd));
126 state_->iterateVec->set(x);
133 state_->stepVec->axpy(-one,
state_->gradientVec->dual());
135 state_->stepVec->axpy(-one,x);
143 if (
state_->searchSize <=
static_cast<Real
>(0) ) {
148 rcon_ = makePtr<ReducedLinearConstraint<Real>>(
proj_->getLinearConstraint(),
151 ns_ = makePtr<NullSpaceOperator<Real>>(
rcon_,x,
152 *
proj_->getResidual());
156template<
typename Real>
161 std::ostream &outStream ) {
164 Real gfnorm(0), gfnormf(0), tol(0), stol(0), snorm(0);
165 Real ftrial(0), pRed(0), rho(1), q(0), delta(0);
166 int flagCG(0), iterCG(0), maxit(0);
169 Ptr<Vector<Real>> s = x.
clone();
170 Ptr<Vector<Real>> gmod = g.
clone(), gfree = g.
clone();
171 Ptr<Vector<Real>> pwa1 = x.
clone(), pwa2 = x.
clone(), pwa3 = x.
clone();
172 Ptr<Vector<Real>> dwa1 = g.
clone(), dwa2 = g.
clone(), dwa3 = g.
clone();
174 Real rhoNM(0), sigmac(0), sigmar(0);
190 *
model_,*dwa1,*dwa2,outStream);
193 delta =
state_->searchSize - snorm;
198 gmod->plus(*
state_->gradientVec);
201 pwa1->set(gfree->dual());
203 gfree->set(pwa1->dual());
206 gfnorm = pwa1->norm();
209 gfnorm = gfree->norm();
213 outStream <<
" Norm of free gradient components: " << gfnorm << std::endl;
214 outStream << std::endl;
219 for (
int i = 0; i <
minit_; ++i) {
221 flagCG = 0; iterCG = 0;
225 if (gfnorm >
zero && delta >
zero) {
226 snorm =
dtrpcg(*s,flagCG,iterCG,*gfree,x,
227 delta,*
model_,bnd,tol,stol,maxit,
228 *pwa1,*dwa1,*pwa2,*dwa2,*pwa3,*dwa3,outStream);
232 outStream <<
" Computation of CG step" << std::endl;
233 outStream <<
" Current face (i): " << i << std::endl;
234 outStream <<
" CG step length: " << snorm << std::endl;
235 outStream <<
" Number of CG iterations: " << iterCG << std::endl;
236 outStream <<
" CG flag: " << flagCG << std::endl;
237 outStream <<
" Total number of iterations: " <<
SPiter_ << std::endl;
238 outStream << std::endl;
246 state_->stepVec->plus(*s);
252 pwa1->set(gfree->dual());
254 gfree->set(pwa1->dual());
257 gfnormf = pwa1->norm();
260 gfnormf = gfree->norm();
263 outStream <<
" Norm of free gradient components: " << gfnormf << std::endl;
264 outStream << std::endl;
269 if (gfnormf <= tol) {
277 else if (flagCG == 2) {
281 else if (delta <=
zero) {
293 ftrial = obj.
value(x,tol0);
302 rho = (rho < rhoNM ? rhoNM : rho );
326 sigmac += pRed; sigmar += pRed;
327 if (ftrial < fmin) { fmin = ftrial; fc = fmin; sigmac =
zero; L = 0; }
330 if (ftrial > fc) { fc = ftrial; sigmac =
zero; }
331 if (L ==
storageNM_) { fr = fc; sigmar = sigmac; }
337 dwa1->set(*
state_->gradientVec);
341 state_->iterateVec->set(x);
353template<
typename Real>
356 std::ostream &outStream)
const {
359 s.
axpy(
static_cast<Real
>(-1),x);
363template<
typename Real>
372 std::ostream &outStream) {
373 const Real half(0.5);
376 Real gs(0), snorm(0);
378 snorm =
dgpstep(s,g,x,-alpha,outStream);
385 q = half * s.
apply(dwa) + gs;
386 interp = (q >
mu0_*gs);
392 while (search && cnt <
redlim_) {
394 snorm =
dgpstep(s,g,x,-alpha,outStream);
398 q = half * s.
apply(dwa) + gs;
399 search = (q >
mu0_*gs);
404 outStream <<
"Cauchy point: The interpolation limit was met without producing sufficient decrease." << std::endl;
405 outStream <<
" Lin-More trust-region algorithm may not converge!" << std::endl;
415 snorm =
dgpstep(s,g,x,-alpha,outStream);
416 if (snorm <= del && cnt <
explim_) {
419 q = half * s.
apply(dwa) + gs;
420 if (q <=
mu0_*gs && std::abs(q-qs) >
qtol_*std::abs(qs)) {
438 snorm =
dgpstep(s,g,x,-alpha,outStream);
441 outStream <<
" Cauchy point" << std::endl;
442 outStream <<
" Step length (alpha): " << alpha << std::endl;
443 outStream <<
" Step length (alpha*g): " << snorm << std::endl;
444 outStream <<
" Model decrease (pRed): " << -q << std::endl;
446 outStream <<
" Number of extrapolation steps: " << cnt << std::endl;
452template<
typename Real>
458 std::ostream &outStream) {
459 const Real
zero(0.0), half(0.5);
461 Real beta(1), snorm(0), gs(0);
467 snorm =
dgpstep(pwa,s,x,beta,outStream);
471 q = half * pwa.
apply(dwa) + gs;
482 outStream << std::endl;
483 outStream <<
" Projected search" << std::endl;
484 outStream <<
" Step length (beta): " << beta << std::endl;
485 outStream <<
" Step length (beta*s): " << snorm << std::endl;
486 outStream <<
" Model decrease (pRed): " << -q << std::endl;
487 outStream <<
" Number of steps: " << nsteps << std::endl;
492template<
typename Real>
496 const Real del)
const {
499 Real rad = ptx*ptx + ptp*(dsq-xtx);
500 rad = std::sqrt(std::max(rad,
zero));
503 sigma = (dsq-xtx)/(ptx+rad);
505 else if (rad >
zero) {
506 sigma = (rad-ptx)/ptp;
514template<
typename Real>
519 const Real tol,
const Real stol,
const int itermax,
522 std::ostream &outStream)
const {
528 const Real
zero(0), one(1), two(2);
529 Real rho(0), kappa(0), beta(0), sigma(0), alpha(0);
530 Real rtr(0), tnorm(0), sMs(0), pMp(0), sMp(0);
544 for (iter = 0; iter < itermax; ++iter) {
550 alpha = (kappa>
zero) ? rho/kappa :
zero;
551 sigma =
dtrqsol(sMs,pMp,sMp,del);
553 if (kappa <= zero || alpha >= sigma) {
556 iflag = (kappa<=
zero) ? 2 : 3;
567 if (rtr <= stol*stol || tnorm <= tol) {
568 sMs = sMs + two*alpha*sMp + alpha*alpha*pMp;
580 sMs = sMs + two*alpha*sMp + alpha*alpha*pMp;
581 sMp = beta*(sMp + alpha*pMp);
582 pMp = (!
hasEcon_ ? rho : p.
dot(p)) + beta*beta*pMp;
585 if (iter == itermax) {
591 return std::sqrt(sMs);
594template<
typename Real>
615template<
typename Real>
638 rcon_->setX(makePtrFromRef(x));
641 ns_->apply(hv,pwa,tol);
645template<
typename Real>
647 std::ios_base::fmtflags osFlags(os.flags());
649 os << std::string(114,
'-') << std::endl;
650 os <<
" Lin-More trust-region method status output definitions" << std::endl << std::endl;
651 os <<
" iter - Number of iterates (steps taken)" << std::endl;
652 os <<
" value - Objective function value" << std::endl;
653 os <<
" gnorm - Norm of the gradient" << std::endl;
654 os <<
" snorm - Norm of the step (update to optimization vector)" << std::endl;
655 os <<
" delta - Trust-Region radius" << std::endl;
656 os <<
" #fval - Number of times the objective function was evaluated" << std::endl;
657 os <<
" #grad - Number of times the gradient was computed" << std::endl;
658 os <<
" #hess - Number of times the Hessian was applied" << std::endl;
659 os <<
" #proj - Number of times the projection was applied" << std::endl;
661 os <<
" tr_flag - Trust-Region flag" << std::endl;
668 os <<
" iterCG - Number of Truncated CG iterations" << std::endl << std::endl;
669 os <<
" flagGC - Trust-Region Truncated CG flag" << std::endl;
675 os << std::string(114,
'-') << std::endl;
678 os << std::setw(6) << std::left <<
"iter";
679 os << std::setw(15) << std::left <<
"value";
680 os << std::setw(15) << std::left <<
"gnorm";
681 os << std::setw(15) << std::left <<
"snorm";
682 os << std::setw(15) << std::left <<
"delta";
683 os << std::setw(10) << std::left <<
"#fval";
684 os << std::setw(10) << std::left <<
"#grad";
685 os << std::setw(10) << std::left <<
"#hess";
686 os << std::setw(10) << std::left <<
"#proj";
687 os << std::setw(10) << std::left <<
"tr_flag";
689 os << std::setw(10) << std::left <<
"iterCG";
690 os << std::setw(10) << std::left <<
"flagCG";
696template<
typename Real>
698 std::ios_base::fmtflags osFlags(os.flags());
699 os << std::endl <<
"Lin-More Trust-Region Method (Type B, Bound Constraints)" << std::endl;
703template<
typename Real>
705 std::ios_base::fmtflags osFlags(os.flags());
706 os << std::scientific << std::setprecision(6);
709 if (
state_->iter == 0 ) {
711 os << std::setw(6) << std::left <<
state_->iter;
712 os << std::setw(15) << std::left <<
state_->value;
713 os << std::setw(15) << std::left <<
state_->gnorm;
714 os << std::setw(15) << std::left <<
"---";
715 os << std::setw(15) << std::left <<
state_->searchSize;
716 os << std::setw(10) << std::left <<
state_->nfval;
717 os << std::setw(10) << std::left <<
state_->ngrad;
718 os << std::setw(10) << std::left <<
nhess_;
719 os << std::setw(10) << std::left <<
state_->nproj;
720 os << std::setw(10) << std::left <<
"---";
722 os << std::setw(10) << std::left <<
"---";
723 os << std::setw(10) << std::left <<
"---";
729 os << std::setw(6) << std::left <<
state_->iter;
730 os << std::setw(15) << std::left <<
state_->value;
731 os << std::setw(15) << std::left <<
state_->gnorm;
732 os << std::setw(15) << std::left <<
state_->snorm;
733 os << std::setw(15) << std::left <<
state_->searchSize;
734 os << std::setw(10) << std::left <<
state_->nfval;
735 os << std::setw(10) << std::left <<
state_->ngrad;
736 os << std::setw(10) << std::left <<
nhess_;
737 os << std::setw(10) << std::left <<
state_->nproj;
738 os << std::setw(10) << std::left <<
TRflag_;
740 os << std::setw(10) << std::left <<
SPiter_;
741 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.
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.
Provides an interface to check status of optimization algorithms.
Provides the interface to evaluate trust-region model functions.
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
virtual void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
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_
bool useSecantHessVec_
Flag to use secant as Hessian (default: false).
Real alpha_
Initial Cauchy point step length (default: 1.0).
Real interpfPS_
Backtracking rate for projected search (default: 0.5).
Ptr< ReducedLinearConstraint< Real > > rcon_
Equality constraint restricted to current active variables.
Real eta1_
Radius decrease threshold (default: 0.05).
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout)
Real TRsafe_
Safeguard size for numerically evaluating ratio (default: 1e2).
void applyFreeHessian(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &pwa) const
Real dcauchy(Vector< Real > &s, Real &alpha, Real &q, const Vector< Real > &x, const Vector< Real > &g, const Real del, TrustRegionModel_U< Real > &model, Vector< Real > &dwa, Vector< Real > &dwa1, std::ostream &outStream=std::cout)
bool writeHeader_
Flag to write header at every iteration.
void writeOutput(std::ostream &os, bool write_header=false) const override
Print iterate status.
bool hasEcon_
Flag signifies if equality constraints exist.
Real tol1_
Absolute tolerance for truncated CG (default: 1e-4).
Real spexp_
Relative tolerance exponent for subproblem solve (default: 1, range: [1,2]).
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 dtrqsol(const Real xtx, const Real ptp, const Real ptx, const Real del) const
int nhess_
Number of Hessian applications.
Real mu0_
Sufficient decrease parameter (default: 1e-2).
Real dtrpcg(Vector< Real > &w, int &iflag, int &iter, const Vector< Real > &g, const Vector< Real > &x, const Real del, TrustRegionModel_U< Real > &model, 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, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &dwa, Vector< Real > &pwa) const
Real eta0_
Step acceptance threshold (default: 0.05).
int SPflag_
Subproblem solver termination flag.
Ptr< TrustRegionModel_U< Real > > model_
Container for trust-region model.
Real interpf_
Backtracking rate for Cauchy point computation (default: 1e-1).
int minit_
Maximum number of minor (subproblem solve) iterations (default: 10).
int maxit_
Maximum number of CG iterations (default: 20).
bool normAlpha_
Normalize initial Cauchy point step length (default: false).
bool interpRad_
Interpolate the trust-region radius if ratio is negative (default: false).
Real qtol_
Relative tolerance for computed decrease in Cauchy point computation (default: 1-8).
Real tol2_
Relative tolerance for truncated CG (default: 1e-2).
Real gamma1_
Radius decrease rate (positive rho) (default: 0.25).
Real dprsrch(Vector< Real > &x, Vector< Real > &s, Real &q, const Vector< Real > &g, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Vector< Real > &pwa, Vector< Real > &dwa, std::ostream &outStream=std::cout)
Real gamma2_
Radius increase rate (default: 2.5).
Real gamma0_
Radius decrease rate (negative rho) (default: 0.0625).
Real eps_
Safeguard for numerically evaluating ratio.
Real delMax_
Maximum trust-region radius (default: ROL_INF).
void writeName(std::ostream &os) const override
Print step name.
unsigned verbosity_
Output level (default: 0).
ESecant esec_
Secant type (default: Limited-Memory BFGS).
Real eta2_
Radius increase threshold (default: 0.9).
Real extrapf_
Extrapolation rate for Cauchy point computation (default: 1e1).
TRUtils::ETRFlag TRflag_
Trust-region exit flag.
bool useSecantPrecond_
Flag to use secant as a preconditioner (default: false).
int redlim_
Maximum number of Cauchy point reduction steps (default: 10).
Real dgpstep(Vector< Real > &s, const Vector< Real > &w, const Vector< Real > &x, const Real alpha, std::ostream &outStream=std::cout) const
int pslim_
Maximum number of projected search steps (default: 20).
Ptr< NullSpaceOperator< Real > > ns_
Null space projection onto reduced equality constraint Jacobian.
LinMoreAlgorithm(ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr)
int SPiter_
Subproblem solver iteration count.
void writeHeader(std::ostream &os) const override
Print iterate header.
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 .
void analyzeRatio(Real &rho, ETRFlag &flag, const Real fold, const Real ftrial, const Real pRed, const Real epsi, std::ostream &outStream=std::cout, const bool print=false)
std::string ETRFlagToString(ETRFlag trf)
Real interpolateRadius(const Vector< Real > &g, const Vector< Real > &s, const Real snorm, const Real pRed, const Real fold, const Real ftrial, const Real del, const Real gamma0, const Real gamma1, const Real eta2, std::ostream &outStream=std::cout, const bool print=false)
std::string NumberToString(T Number)
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
ESecant StringToESecant(std::string s)
Real ROL_OVERFLOW(void)
Platform-dependent maximum double.
std::string ECGFlagToString(ECGFlag cgf)