10#ifndef ROL_TYPEB_COLEMANLIALGORITHM_DEF_HPP
11#define ROL_TYPEB_COLEMANLIALGORITHM_DEF_HPP
16template<
typename Real>
23 ParameterList &trlist = list.sublist(
"Step").sublist(
"Trust Region");
25 state_->searchSize = trlist.get(
"Initial Radius", -1.0);
27 eta0_ = trlist.get(
"Step Acceptance Threshold", 0.05);
28 eta1_ = trlist.get(
"Radius Shrinking Threshold", 0.05);
29 eta2_ = trlist.get(
"Radius Growing Threshold", 0.9);
30 gamma0_ = trlist.get(
"Radius Shrinking Rate (Negative rho)", 0.0625);
31 gamma1_ = trlist.get(
"Radius Shrinking Rate (Positive rho)", 0.25);
32 gamma2_ = trlist.get(
"Radius Growing Rate", 2.5);
33 TRsafe_ = trlist.get(
"Safeguard Size", 100.0);
35 interpRad_ = trlist.get(
"Use Radius Interpolation",
false);
37 storageNM_ = trlist.get(
"Nonmonotone Storage Size", 0);
40 maxit_ = list.sublist(
"General").sublist(
"Krylov").get(
"Iteration Limit", 20);
41 tol1_ = list.sublist(
"General").sublist(
"Krylov").get(
"Absolute Tolerance", 1e-4);
42 tol2_ = list.sublist(
"General").sublist(
"Krylov").get(
"Relative Tolerance", 1e-2);
44 ROL::ParameterList &lmlist = trlist.sublist(
"Coleman-Li");
45 mu0_ = lmlist.get(
"Sufficient Decrease Parameter", 1e-2);
46 spexp_ = lmlist.get(
"Relative Tolerance Exponent", 1.0);
47 spexp_ = std::max(
static_cast<Real
>(1),std::min(
spexp_,
static_cast<Real
>(2)));
48 alphaMax_ = lmlist.get(
"Relaxation Safeguard", 0.999);
51 verbosity_ = list.sublist(
"General").get(
"Output Level",0);
54 useSecantPrecond_ = list.sublist(
"General").sublist(
"Secant").get(
"Use as Preconditioner",
false);
55 useSecantHessVec_ = list.sublist(
"General").sublist(
"Secant").get(
"Use as Hessian",
false);
60 model_ = makePtr<TrustRegionModel_U<Real>>(list,secant,mode);
61 if (secant == nullPtr) {
62 std::string secantType = list.sublist(
"General").sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS");
67template<
typename Real>
72 std::ostream &outStream) {
75 if (
proj_ == nullPtr) {
76 proj_ = makePtr<PolyhedralProjection<Real>>(makePtrFromRef(bnd));
84 proj_->getBoundConstraint()->projectInterior(x);
state_->nproj++;
85 state_->iterateVec->set(x);
94 state_->stepVec->axpy(-one,x);
98 if (
state_->searchSize <=
static_cast<Real
>(0) ) {
103 rcon_ = makePtr<ReducedLinearConstraint<Real>>(
proj_->getLinearConstraint(),
106 ns_ = makePtr<NullSpaceOperator<Real>>(
rcon_,x,
107 *
proj_->getResidual());
111template<
typename Real>
116 std::ostream &outStream ) {
117 const Real
zero(0), one(1), half(0.5);
119 Real tol(0), stol(0), snorm(0);
120 Real ftrial(0), pRed(0), rho(1), alpha(1);
123 Ptr<Vector<Real>> pwa1 = x.
clone(), pwa2 = x.
clone(), pwa3 = x.
clone();
124 Ptr<Vector<Real>> pwa4 = x.
clone(), pwa5 = x.
clone();
125 Ptr<Vector<Real>> dwa1 = g.
clone(), dwa2 = g.
clone(), dwa3 = g.
clone();
127 Real rhoNM(0), sigmac(0), sigmar(0), sBs(0), gs(0);
146 pwa5->set(
state_->gradientVec->dual());
149 *pwa1,*dwa1,*pwa2,*dwa2,*pwa3,*pwa4,*dwa3,outStream);
151 outStream <<
" Computation of CG step" << std::endl;
152 outStream <<
" CG step length: " << snorm << std::endl;
153 outStream <<
" Number of CG iterations: " <<
SPiter_ << std::endl;
154 outStream <<
" CG flag: " <<
SPflag_ << std::endl;
155 outStream << std::endl;
162 state_->stepVec->set(*pwa1);
163 state_->snorm = alpha * snorm;
169 sBs = dwa1->apply(*
state_->stepVec);
170 pRed = - half * sBs - gs;
174 ftrial = obj.
value(x,tol0);
183 rho = (rho < rhoNM ? rhoNM : rho );
207 sigmac += pRed; sigmar += pRed;
208 if (ftrial < fmin) { fmin = ftrial; fc = fmin; sigmac =
zero; L = 0; }
211 if (ftrial > fc) { fc = ftrial; sigmac =
zero; }
212 if (L ==
storageNM_) { fr = fc; sigmar = sigmac; }
218 dwa1->set(*
state_->gradientVec);
222 state_->iterateVec->set(x);
234template<
typename Real>
237 std::ostream &outStream)
const {
239 proj_->getBoundConstraint()->projectInterior(s);
state_->nproj++;
240 s.
axpy(
static_cast<Real
>(-1),x);
244template<
typename Real>
248 const Real del)
const {
251 Real rad = ptx*ptx + ptp*(dsq-xtx);
252 rad = std::sqrt(std::max(rad,
zero));
255 sigma = (dsq-xtx)/(ptx+rad);
257 else if (rad >
zero) {
258 sigma = (rad-ptx)/ptp;
266template<
typename Real>
272 const Real tol,
const Real stol,
276 std::ostream &outStream)
const {
282 const Real
zero(0), one(1), two(2);
283 Real rho(0), kappa(0), beta(0), sigma(0), alpha(0);
284 Real rtr(0), tnorm(0), sMs(0), pMp(0), sMp(0);
298 for (iter = 0; iter <
maxit_; ++iter) {
304 alpha = (kappa>
zero) ? rho/kappa :
zero;
305 sigma =
dtrqsol(sMs,pMp,sMp,del);
307 if (kappa <= zero || alpha >= sigma) {
310 iflag = (kappa<=
zero) ? 2 : 3;
321 if (rtr <= stol*stol || tnorm <= tol) {
322 sMs = sMs + two*alpha*sMp + alpha*alpha*pMp;
334 sMs = sMs + two*alpha*sMp + alpha*alpha*pMp;
335 sMp = beta*(sMp + alpha*pMp);
336 pMp = (!
hasEcon_ ? rho : p.
dot(p)) + beta*beta*pMp;
345 return std::sqrt(sMs);
348template<
typename Real>
359template<
typename Real>
370 applyC(pwa2,v,x,g,bnd,pwa1);
374template<
typename Real>
387template<
typename Real>
389 std::ios_base::fmtflags osFlags(os.flags());
391 os << std::string(114,
'-') << std::endl;
392 os <<
" Coleman-Li affine-scaling trust-region method status output definitions" << std::endl << std::endl;
393 os <<
" iter - Number of iterates (steps taken)" << std::endl;
394 os <<
" value - Objective function value" << std::endl;
395 os <<
" gnorm - Norm of the gradient" << std::endl;
396 os <<
" snorm - Norm of the step (update to optimization vector)" << std::endl;
397 os <<
" delta - Trust-Region radius" << std::endl;
398 os <<
" #fval - Number of times the objective function was evaluated" << std::endl;
399 os <<
" #grad - Number of times the gradient was computed" << std::endl;
400 os <<
" #hess - Number of times the Hessian was applied" << std::endl;
401 os <<
" #proj - Number of times the projection was applied" << std::endl;
403 os <<
" tr_flag - Trust-Region flag" << std::endl;
409 os <<
" iterCG - Number of Truncated CG iterations" << std::endl << std::endl;
410 os <<
" flagGC - Trust-Region Truncated CG flag" << std::endl;
415 os << std::string(114,
'-') << std::endl;
418 os << std::setw(6) << std::left <<
"iter";
419 os << std::setw(15) << std::left <<
"value";
420 os << std::setw(15) << std::left <<
"gnorm";
421 os << std::setw(15) << std::left <<
"snorm";
422 os << std::setw(15) << std::left <<
"delta";
423 os << std::setw(10) << std::left <<
"#fval";
424 os << std::setw(10) << std::left <<
"#grad";
425 os << std::setw(10) << std::left <<
"#hess";
426 os << std::setw(10) << std::left <<
"#proj";
427 os << std::setw(10) << std::left <<
"tr_flag";
428 os << std::setw(10) << std::left <<
"iterCG";
429 os << std::setw(10) << std::left <<
"flagCG";
434template<
typename Real>
436 std::ios_base::fmtflags osFlags(os.flags());
437 os << std::endl <<
"Coleman-Li Affine-Scaling Trust-Region Method (Type B, Bound Constraints)" << std::endl;
441template<
typename Real>
443 std::ios_base::fmtflags osFlags(os.flags());
444 os << std::scientific << std::setprecision(6);
447 if (
state_->iter == 0 ) {
449 os << std::setw(6) << std::left <<
state_->iter;
450 os << std::setw(15) << std::left <<
state_->value;
451 os << std::setw(15) << std::left <<
state_->gnorm;
452 os << std::setw(15) << std::left <<
"---";
453 os << std::setw(15) << std::left <<
state_->searchSize;
454 os << std::setw(10) << std::left <<
state_->nfval;
455 os << std::setw(10) << std::left <<
state_->ngrad;
456 os << std::setw(10) << std::left <<
nhess_;
457 os << std::setw(10) << std::left <<
state_->nproj;
458 os << std::setw(10) << std::left <<
"---";
459 os << std::setw(10) << std::left <<
"---";
460 os << std::setw(10) << std::left <<
"---";
465 os << std::setw(6) << std::left <<
state_->iter;
466 os << std::setw(15) << std::left <<
state_->value;
467 os << std::setw(15) << std::left <<
state_->gnorm;
468 os << std::setw(15) << std::left <<
state_->snorm;
469 os << std::setw(15) << std::left <<
state_->searchSize;
470 os << std::setw(10) << std::left <<
state_->nfval;
471 os << std::setw(10) << std::left <<
state_->ngrad;
472 os << std::setw(10) << std::left <<
nhess_;
473 os << std::setw(10) << std::left <<
state_->nproj;
474 os << std::setw(10) << std::left <<
TRflag_;
475 os << std::setw(10) << std::left <<
SPiter_;
476 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 void applyInverseScalingFunction(Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const
Apply inverse scaling function.
virtual void applyScalingFunctionJacobian(Vector< Real > &dv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g) const
Apply scaling function Jacobian.
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_
Real dgpstep(Vector< Real > &s, const Vector< Real > &w, const Vector< Real > &x, const Real alpha, std::ostream &outStream=std::cout) const
int nhess_
Number of Hessian applications.
bool useSecantPrecond_
Flag to use secant as a preconditioner (default: false).
int maxit_
Maximum number of CG iterations (default: 20).
void writeOutput(std::ostream &os, const bool write_header=false) const override
Print iterate status.
Ptr< ReducedLinearConstraint< Real > > rcon_
Equality constraint restricted to current active variables.
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout)
bool writeHeader_
Flag to write header at every iteration.
Real dtrpcg(Vector< Real > &w, int &iflag, int &iter, const Vector< Real > &g, const Vector< Real > &x, const Vector< Real > &gdual, const Real del, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, const Real tol, const Real stol, Vector< Real > &p, Vector< Real > &q, Vector< Real > &r, Vector< Real > &t, Vector< Real > &pwa1, Vector< Real > &pwa2, Vector< Real > &dwa, std::ostream &outStream=std::cout) const
Real dtrqsol(const Real xtx, const Real ptp, const Real ptx, const Real del) const
Real alphaMax_
Maximum value of relaxation parameter (default: 0.999).
void applyPrecond(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &dwa, Vector< Real > &pwa) const
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...
int SPflag_
Subproblem solver termination flag.
ColemanLiAlgorithm(ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr)
int SPiter_
Subproblem solver iteration count.
TRUtils::ETRFlag TRflag_
Trust-region exit flag.
ESecant esec_
Secant type (default: Limited-Memory BFGS).
bool useSecantHessVec_
Flag to use secant as Hessian (default: false).
unsigned verbosity_
Output level (default: 0).
bool interpRad_
Interpolate the trust-region radius if ratio is negative (default: false).
Real spexp_
Relative tolerance exponent for subproblem solve (default: 1, range: [1,2]).
Real gamma0_
Radius decrease rate (negative rho) (default: 0.0625).
Real tol1_
Absolute tolerance for truncated CG (default: 1e-4).
Real eta1_
Radius decrease threshold (default: 0.05).
Real gamma1_
Radius decrease rate (positive rho) (default: 0.25).
Real delMax_
Maximum trust-region radius (default: ROL_INF).
Ptr< NullSpaceOperator< Real > > ns_
Null space projection onto reduced equality constraint Jacobian.
Real eta2_
Radius increase threshold (default: 0.9).
void writeName(std::ostream &os) const override
Print step name.
bool hasEcon_
Flag signifies if equality constraints exist.
Real eta0_
Step acceptance threshold (default: 0.05).
Real TRsafe_
Safeguard size for numerically evaluating ratio (default: 1e2).
Real mu0_
Sufficient decrease parameter (default: 1e-2).
Real gamma2_
Radius increase rate (default: 2.5).
void applyHessian(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, TrustRegionModel_U< Real > &model, BoundConstraint< Real > &bnd, Real &tol, Vector< Real > &pwa1, Vector< Real > &pwa2) const
Real eps_
Safeguard for numerically evaluating ratio.
Ptr< TrustRegionModel_U< Real > > model_
Container for trust-region model.
Real tol2_
Relative tolerance for truncated CG (default: 1e-2).
void applyC(Vector< Real > &Cv, const Vector< Real > &v, const Vector< Real > &x, const Vector< Real > &g, BoundConstraint< Real > &bnd, Vector< Real > &pwa) const
void writeHeader(std::ostream &os) const override
Print iterate header.
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)