44 #ifndef ROL_TYPEB_TRUSTREGIONSPGALGORITHM_DEF_HPP
45 #define ROL_TYPEB_TRUSTREGIONSPGALGORITHM_DEF_HPP
52 template<
typename Real>
59 ParameterList &trlist = list.sublist(
"Step").sublist(
"Trust Region");
61 state_->searchSize = trlist.get(
"Initial Radius", -1.0);
62 delMax_ = trlist.get(
"Maximum Radius", ROL_INF<Real>());
63 eta0_ = trlist.get(
"Step Acceptance Threshold", 0.05);
64 eta1_ = trlist.get(
"Radius Shrinking Threshold", 0.05);
65 eta2_ = trlist.get(
"Radius Growing Threshold", 0.9);
66 gamma0_ = trlist.get(
"Radius Shrinking Rate (Negative rho)", 0.0625);
67 gamma1_ = trlist.get(
"Radius Shrinking Rate (Positive rho)", 0.25);
68 gamma2_ = trlist.get(
"Radius Growing Rate", 2.5);
69 TRsafe_ = trlist.get(
"Safeguard Size", 100.0);
70 eps_ = TRsafe_*ROL_EPSILON<Real>();
71 interpRad_ = trlist.get(
"Use Radius Interpolation",
false);
72 verbosity_ = trlist.sublist(
"General").get(
"Output Level", 0);
74 storageNM_ = trlist.get(
"Nonmonotone Storage Size", 0);
75 useNM_ = (storageNM_ <= 0 ?
false :
true);
77 ROL::ParameterList &lmlist = trlist.sublist(
"SPG");
78 mu0_ = lmlist.get(
"Sufficient Decrease Parameter", 1e-2);
79 spexp_ = lmlist.get(
"Relative Tolerance Exponent", 1.0);
80 spexp_ = std::max(static_cast<Real>(1),std::min(spexp_,static_cast<Real>(2)));
81 redlim_ = lmlist.sublist(
"Cauchy Point").get(
"Maximum Number of Reduction Steps", 10);
82 explim_ = lmlist.sublist(
"Cauchy Point").get(
"Maximum Number of Expansion Steps", 10);
83 alpha_ = lmlist.sublist(
"Cauchy Point").get(
"Initial Step Size", 1.0);
84 normAlpha_ = lmlist.sublist(
"Cauchy Point").get(
"Normalize Initial Step Size",
false);
85 interpf_ = lmlist.sublist(
"Cauchy Point").get(
"Reduction Rate", 0.1);
86 extrapf_ = lmlist.sublist(
"Cauchy Point").get(
"Expansion Rate", 10.0);
87 qtol_ = lmlist.sublist(
"Cauchy Point").get(
"Decrease Tolerance", 1e-8);
89 lambdaMin_ = lmlist.sublist(
"Solver").get(
"Minimum Spectral Step Size", 1e-8);
90 lambdaMax_ = lmlist.sublist(
"Solver").get(
"Maximum Spectral Step Size", 1e8);
91 gamma_ = lmlist.sublist(
"Solver").get(
"Sufficient Decrease Tolerance", 1e-4);
92 maxSize_ = lmlist.sublist(
"Solver").get(
"Maximum Storage Size", 10);
93 maxit_ = lmlist.sublist(
"Solver").get(
"Iteration Limit", 25);
94 tol1_ = lmlist.sublist(
"Solver").get(
"Absolute Tolerance", 1e-4);
95 tol2_ = lmlist.sublist(
"Solver").get(
"Relative Tolerance", 1e-2);
96 useMin_ = lmlist.sublist(
"Solver").get(
"Use Smallest Model Iterate",
true);
97 useNMSP_ = lmlist.sublist(
"Solver").get(
"Use Nonmonotone Search",
false);
98 useSimpleSPG_ = !lmlist.sublist(
"Solver").get(
"Compute Cauchy Point",
true);
100 ParameterList &glist = list.sublist(
"General");
102 useInexact_.push_back(glist.get(
"Inexact Objective Function",
false));
103 useInexact_.push_back(glist.get(
"Inexact Gradient",
false));
104 useInexact_.push_back(glist.get(
"Inexact Hessian-Times-A-Vector",
false));
106 ParameterList &ilist = trlist.sublist(
"Inexact").sublist(
"Gradient");
107 scale0_ = ilist.get(
"Tolerance Scaling", static_cast<Real>(0.1));
108 scale1_ = ilist.get(
"Relative Tolerance", static_cast<Real>(2));
110 ParameterList &vlist = trlist.sublist(
"Inexact").sublist(
"Value");
111 scale_ = vlist.get(
"Tolerance Scaling", static_cast<Real>(1.e-1));
112 omega_ = vlist.get(
"Exponent", static_cast<Real>(0.9));
113 force_ = vlist.get(
"Forcing Sequence Initial Value", static_cast<Real>(1.0));
114 updateIter_ = vlist.get(
"Forcing Sequence Update Frequency", static_cast<int>(10));
115 forceFactor_ = vlist.get(
"Forcing Sequence Reduction Factor", static_cast<Real>(0.1));
117 verbosity_ = list.sublist(
"General").get(
"Output Level",0);
118 writeHeader_ = verbosity_ > 2;
120 useSecantPrecond_ = list.sublist(
"General").sublist(
"Secant").get(
"Use as Preconditioner",
false);
121 useSecantHessVec_ = list.sublist(
"General").sublist(
"Secant").get(
"Use as Hessian",
false);
126 model_ = makePtr<TrustRegionModel_U<Real>>(list,secant,mode);
127 if (secant == nullPtr) {
128 esec_ =
StringToESecant(list.sublist(
"General").sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS"));
132 template<
typename Real>
138 std::ostream &outStream) {
140 if (proj_ == nullPtr)
141 proj_ = makePtr<PolyhedralProjection<Real>>(makePtrFromRef(bnd));
146 proj_->project(x,outStream); state_->nproj++;
147 state_->iterateVec->set(x);
149 state_->value = obj.
value(x,ftol);
152 computeGradient(x,*state_->gradientVec,*state_->stepVec,state_->searchSize,obj,
true,gtol_,state_->gnorm,outStream);
159 state_->snorm = ROL_INF<Real>();
161 if (normAlpha_) alpha_ /= state_->gradientVec->norm();
163 if ( state_->searchSize <= static_cast<Real>(0) )
164 state_->searchSize = state_->gradientVec->norm();
167 template<
typename Real>
176 outTol = std::sqrt(ROL_EPSILON<Real>());
177 if ( useInexact_[0] ) {
178 if (!(iter%updateIter_) && (iter!=0)) force_ *= forceFactor_;
180 Real eta =
static_cast<Real
>(0.999)*std::min(eta1_,one-eta2_);
181 outTol = scale_*std::pow(eta*std::min(pRed,force_),one/omega_);
182 if (inTol > outTol) fold = obj.
value(xold,outTol);
186 Real fval = obj.
value(x,outTol);
190 template<
typename Real>
199 std::ostream &outStream)
const {
200 if ( useInexact_[1] ) {
202 Real gtol0 = scale0_*del;
203 if (accept) gtol = gtol0 + one;
204 else gtol0 = scale0_*std::min(gnorm,del);
205 while ( gtol > gtol0 ) {
209 gtol0 = scale0_*std::min(gnorm,del);
214 gtol = std::sqrt(ROL_EPSILON<Real>());
221 template<
typename Real>
226 std::ostream &outStream ) {
227 const Real
zero(0), one(1);
229 Real inTol =
static_cast<Real
>(0.1)*ROL_OVERFLOW<Real>(), outTol(inTol);
230 Real ftrial(0), pRed(0), rho(1), q(0);
232 std::vector<std::string> output;
233 initialize(x,g,inTol,obj,bnd,outStream);
234 Ptr<Vector<Real>> gmod = g.
clone();
235 Ptr<Vector<Real>> pwa1 = x.
clone(), pwa2 = x.
clone();
236 Ptr<Vector<Real>> pwa3 = x.
clone(), pwa4 = x.
clone();
237 Ptr<Vector<Real>> pwa5 = x.
clone(), pwa6 = x.
clone();
238 Ptr<Vector<Real>> pwa7 = x.
clone();
239 Ptr<Vector<Real>> dwa1 = g.
clone(), dwa2 = g.
clone();
241 Real rhoNM(0), sigmac(0), sigmar(0);
242 Real fr(state_->value), fc(state_->value), fmin(state_->value);
247 if (verbosity_ > 0) writeOutput(outStream,
true);
249 while (status_->check(*state_)) {
251 model_->setData(obj,*state_->iterateVec,*state_->gradientVec,gtol_);
255 gmod->set(*state_->gradientVec);
257 dpsg_simple(x,q,*gmod,*state_->iterateVec,state_->searchSize,*model_,
258 *pwa1,*pwa2,*dwa1,outStream);
261 dcauchy(*state_->stepVec,alpha_,q,*state_->iterateVec,
262 state_->gradientVec->dual(),state_->searchSize,
263 *model_,*dwa1,*dwa2,outStream);
264 x.
plus(*state_->stepVec);
270 dpsg(x,q,*gmod,*state_->iterateVec,state_->searchSize,*model_,
271 *pwa1,*pwa2,*pwa3,*pwa4,*pwa5,*pwa6,*pwa7,*dwa1,outStream);
276 state_->stepVec->set(x); state_->stepVec->axpy(-one,*state_->iterateVec);
277 state_->snorm = state_->stepVec->norm();
280 ftrial = computeValue(inTol,outTol,pRed,state_->value,state_->iter,x,*state_->iterateVec,obj);
285 TRUtils::analyzeRatio<Real>(rho,TRflag_,state_->value,ftrial,pRed,eps_,outStream,verbosity_>1);
287 TRUtils::analyzeRatio<Real>(rhoNM,TRflagNM,fr,ftrial,pRed+sigmar,eps_,outStream,verbosity_>1);
288 TRflag_ = (rho < rhoNM ? TRflagNM : TRflag_);
289 rho = (rho < rhoNM ? rhoNM : rho );
296 x.
set(*state_->iterateVec);
300 state_->searchSize = TRUtils::interpolateRadius<Real>(*state_->gradientVec,*state_->stepVec,
301 state_->snorm,pRed,state_->value,ftrial,state_->searchSize,gamma0_,gamma1_,eta2_,
302 outStream,verbosity_>1);
305 state_->searchSize = gamma1_*std::min(state_->snorm,state_->searchSize);
307 computeGradient(x,*state_->gradientVec,*pwa1,state_->searchSize,obj,
false,gtol_,state_->gnorm,outStream);
311 state_->value = ftrial;
315 sigmac += pRed; sigmar += pRed;
316 if (ftrial < fmin) { fmin = ftrial; fc = fmin; sigmac =
zero; L = 0; }
319 if (ftrial > fc) { fc = ftrial; sigmac =
zero; }
320 if (L == storageNM_) { fr = fc; sigmar = sigmac; }
324 if (rho >= eta2_) state_->searchSize = std::min(gamma2_*state_->searchSize, delMax_);
326 dwa1->set(*state_->gradientVec);
327 computeGradient(x,*state_->gradientVec,*pwa1,state_->searchSize,obj,
true,gtol_,state_->gnorm,outStream);
329 state_->iterateVec->set(x);
331 model_->update(x,*state_->stepVec,*dwa1,*state_->gradientVec,
332 state_->snorm,state_->iter);
336 if (verbosity_ > 0) writeOutput(outStream,writeHeader_);
341 template<
typename Real>
344 std::ostream &outStream)
const {
346 proj_->project(s,outStream); state_->nproj++;
347 s.
axpy(static_cast<Real>(-1),x);
351 template<
typename Real>
360 std::ostream &outStream) {
361 const Real half(0.5);
363 Real tol = std::sqrt(ROL_EPSILON<Real>());
365 Real gs(0), snorm(0);
367 snorm = dgpstep(s,g,x,-alpha,outStream);
372 model.
hessVec(dwa,s,x,tol); nhess_++;
375 q = half * s.
apply(dwa) + gs;
376 interp = (q > mu0_*gs);
384 snorm = dgpstep(s,g,x,-alpha,outStream);
386 model.
hessVec(dwa,s,x,tol); nhess_++;
389 q = half * s.
apply(dwa) + gs;
390 search = (q > mu0_*gs) && (cnt < redlim_);
402 snorm = dgpstep(s,g,x,-alpha,outStream);
403 if (snorm <= del && cnt < explim_) {
404 model.
hessVec(dwa,s,x,tol); nhess_++;
407 q = half * s.
apply(dwa) + gs;
408 if (q <= mu0_*gs && std::abs(q-qs) > qtol_*std::abs(qs)) {
426 snorm = dgpstep(s,g,x,-alpha,outStream);
428 if (verbosity_ > 1) {
429 outStream <<
" Cauchy point" << std::endl;
430 outStream <<
" Step length (alpha): " << alpha << std::endl;
431 outStream <<
" Step length (alpha*g): " << snorm << std::endl;
432 outStream <<
" Model decrease (pRed): " << -q << std::endl;
434 outStream <<
" Number of extrapolation steps: " << cnt << std::endl;
440 template<
typename Real>
450 std::ostream &outStream) {
458 const Real half(0.5), one(1), safeguard(1e2*ROL_EPSILON<Real>());
459 Real tol(std::sqrt(ROL_EPSILON<Real>()));
460 Real alpha(1), alphaMax(1), s0s0(0), ss0(0), sHs(0), lambdaTmp(1), snorm(0);
467 Real coeff = one/gmod.
norm();
468 Real lambda = std::max(lambdaMin_,std::min(coeff,lambdaMax_));
470 proj_->project(pwa,outStream); state_->nproj++;
472 Real gs = gmod.
apply(pwa);
473 Real ss = pwa.
dot(pwa);
474 Real gnorm = std::sqrt(ss);
477 const Real gtol = std::min(tol1_,tol2_*gnorm);
480 outStream <<
" Spectral Projected Gradient" << std::endl;
483 while (SPiter_ < maxit_) {
487 model.
hessVec(dwa,pwa,x,tol); nhess_++;
488 sHs = dwa.
apply(pwa);
492 if (gnorm >= del-safeguard) {
494 alphaMax = std::min(one, (-ss0 + std::sqrt(ss0*ss0 - ss*(s0s0-del*del)))/ss);
496 if (sHs <= safeguard)
499 alpha = std::min(alphaMax, -gs/sHs);
502 q += alpha * (gs + half * alpha * sHs);
503 gmod.
axpy(alpha,dwa);
507 pwa1.
set(y); pwa1.
axpy(-one,x);
508 s0s0 = pwa1.
dot(pwa1);
509 snorm = std::sqrt(s0s0);
511 if (verbosity_ > 1) {
512 outStream << std::endl;
513 outStream <<
" Iterate: " << SPiter_ << std::endl;
514 outStream <<
" Spectral step length (lambda): " << lambda << std::endl;
515 outStream <<
" Step length (alpha): " << alpha << std::endl;
516 outStream <<
" Model decrease (pRed): " << -q << std::endl;
517 outStream <<
" Optimality criterion: " << gnorm << std::endl;
518 outStream <<
" Step norm: " << snorm << std::endl;
519 outStream << std::endl;
522 if (snorm >= del - safeguard) { SPflag_ = 2;
break; }
525 lambdaTmp = (sHs <= safeguard ? one/gmod.
norm() : ss/sHs);
526 lambda = std::max(lambdaMin_,std::min(lambdaTmp,lambdaMax_));
528 proj_->project(pwa,outStream); state_->nproj++;
530 gs = gmod.
apply(pwa);
532 gnorm = std::sqrt(ss);
534 if (gnorm <= gtol) { SPflag_ = 0;
break; }
536 SPflag_ = (SPiter_==maxit_) ? 1 : SPflag_;
539 template<
typename Real>
554 std::ostream &outStream) {
562 const Real
zero(0), half(0.5), one(1), two(2);
563 Real tol(std::sqrt(ROL_EPSILON<Real>()));
564 Real alpha(1), sHs(0), alphaTmp(1), mmax(0), qmin(0), lambdaTmp(1);
565 std::deque<Real> mqueue; mqueue.push_back(q);
567 if (useNMSP_ && useMin_) { qmin = q; ymin.
set(y); }
571 pwa.
set(y); pwa.
axpy(-one,pwa1);
572 dproj(pwa,x,del,pwa2,pwa3,pwa4,pwa5,outStream);
574 Real gnorm = pwa.
norm();
575 const Real gtol = std::min(tol1_,tol2_*gnorm);
578 Real coeff = one/gmod.
norm();
579 Real lambda = std::max(lambdaMin_,std::min(coeff,lambdaMax_));
580 pwa.
set(y); pwa.
axpy(-lambda,pwa1);
581 dproj(pwa,x,del,pwa2,pwa3,pwa4,pwa5,outStream);
583 Real gs = gmod.
apply(pwa);
584 Real ss = pwa.
dot(pwa);
587 outStream <<
" Spectral Projected Gradient" << std::endl;
590 while (SPiter_ < maxit_) {
594 model.
hessVec(dwa,pwa,x,tol); nhess_++;
595 sHs = dwa.
apply(pwa);
599 mmax = *std::max_element(mqueue.begin(),mqueue.end());
600 alphaTmp = (-(one-gamma_)*gs + std::sqrt(std::pow((one-gamma_)*gs,two)-two*sHs*(q-mmax)))/sHs;
605 alpha = (sHs >
zero ? std::min(one,std::max(
zero,alphaTmp)) : one);
608 q += alpha * (gs + half * alpha * sHs);
609 gmod.
axpy(alpha,dwa);
614 if (static_cast<int>(mqueue.size())==maxSize_) mqueue.pop_front();
616 if (useMin_ && q <= qmin) { qmin = q; ymin.
set(y); }
621 pwa.
set(y); pwa.
axpy(-one,pwa1);
622 dproj(pwa,x,del,pwa2,pwa3,pwa4,pwa5,outStream);
626 if (verbosity_ > 1) {
627 outStream << std::endl;
628 outStream <<
" Iterate: " << SPiter_ << std::endl;
629 outStream <<
" Spectral step length (lambda): " << lambda << std::endl;
630 outStream <<
" Step length (alpha): " << alpha << std::endl;
631 outStream <<
" Model decrease (pRed): " << -q << std::endl;
632 outStream <<
" Optimality criterion: " << gnorm << std::endl;
633 outStream << std::endl;
635 if (gnorm < gtol)
break;
639 lambdaTmp = (sHs == 0 ? coeff : ss/sHs);
640 lambda = std::max(lambdaMin_,std::min(lambdaTmp,lambdaMax_));
641 pwa.
set(y); pwa.
axpy(-lambda,pwa1);
642 dproj(pwa,x,del,pwa2,pwa3,pwa4,pwa5,outStream);
644 gs = gmod.
apply(pwa);
647 if (useNMSP_ && useMin_) { q = qmin; y.
set(ymin); }
648 SPflag_ = (SPiter_==maxit_) ? 1 : 0;
651 template<
typename Real>
659 std::ostream &outStream)
const {
661 const Real
zero(0), half(0.5), one(1), two(2), three(3);
662 const Real eps(ROL_EPSILON<Real>()), tol0(1e1*eps), fudge(1.0-1e-2*sqrt(eps));
663 Real f0(0), f1(0), fc(0), t0(0), t1(1), tc(0), d1(1), d2(1), tol(1);
664 Real p(0), q(0), r(0), s(0), m(0);
665 int cnt(state_->nproj);
667 proj_->project(y1,outStream); state_->nproj++;
668 pwa.
set(y1); pwa.
axpy(-one,x0);
675 tc = t0; fc = f0; yc.
set(y0);
679 if (std::abs(fc-del) < std::abs(f1-del)) {
680 t0 = t1; t1 = tc; tc = t0;
681 f0 = f1; f1 = fc; fc = f0;
684 tol = two*eps*std::abs(t1) + half*tol0;
686 if (std::abs(m) <= tol) { code = 1;
break; }
687 if ((f1 >= fudge*del && f1 <= del))
break;
688 if (std::abs(d1) < tol || std::abs(f0-del) <= std::abs(f1-del)) {
692 s = (f1-del)/(f0-del);
698 q = (f0-del)/(fc-del);
699 r = (f1-del)/(fc-del);
700 p = s*(two*m*q*(q-r)-(t1-t0)*(r-one));
701 q = (q-one)*(r-one)*(s-one);
703 if (p >
zero) q = -q;
707 if (two*p < three*m*q-std::abs(tol*q) && p < std::abs(half*s*q)) {
714 t0 = t1; f0 = f1; y0.
set(y1);
715 if (std::abs(d2) > tol) t1 += d2;
716 else if (m >
zero) t1 += tol;
719 proj_->project(y1,outStream); state_->nproj++;
720 pwa.
set(y1); pwa.
axpy(-one,x0);
722 if ((f1 > del && fc > del) || (f1 <= del && fc <= del)) {
723 tc = t0; fc = f0; yc.
set(y0);
727 if (code==1 && f1>del) x.
set(yc);
729 if (verbosity_ > 1) {
730 outStream << std::endl;
731 outStream <<
" Trust-Region Subproblem Projection" << std::endl;
732 outStream <<
" Number of polyhedral projections: " << state_->nproj-cnt << std::endl;
733 if (code == 1 && f1 > del) {
734 outStream <<
" Transformed Multiplier: " << tc << std::endl;
735 outStream <<
" Dual Residual: " << fc-del << std::endl;
738 outStream <<
" Transformed Multiplier: " << t1 << std::endl;
739 outStream <<
" Dual Residual: " << f1-del << std::endl;
741 outStream <<
" Exit Code: " << code << std::endl;
742 outStream << std::endl;
920 template<
typename Real>
922 std::ios_base::fmtflags osFlags(os.flags());
923 if (verbosity_ > 1) {
924 os << std::string(114,
'-') << std::endl;
925 os <<
" SPG trust-region method status output definitions" << std::endl << std::endl;
926 os <<
" iter - Number of iterates (steps taken)" << std::endl;
927 os <<
" value - Objective function value" << std::endl;
928 os <<
" gnorm - Norm of the gradient" << std::endl;
929 os <<
" snorm - Norm of the step (update to optimization vector)" << std::endl;
930 os <<
" delta - Trust-Region radius" << std::endl;
931 os <<
" #fval - Number of times the objective function was evaluated" << std::endl;
932 os <<
" #grad - Number of times the gradient was computed" << std::endl;
933 os <<
" #hess - Number of times the Hessian was applied" << std::endl;
934 os <<
" #proj - Number of times the projection was computed" << std::endl;
936 os <<
" tr_flag - Trust-Region flag" << std::endl;
942 os <<
" iterSPG - Number of Spectral Projected Gradient iterations" << std::endl << std::endl;
943 os <<
" flagSPG - Trust-Region Truncated CG flag" << std::endl;
944 os <<
" 0 - Converged" << std::endl;
945 os <<
" 1 - Iteration Limit Exceeded" << std::endl;
946 os << std::string(114,
'-') << std::endl;
949 os << std::setw(6) << std::left <<
"iter";
950 os << std::setw(15) << std::left <<
"value";
951 os << std::setw(15) << std::left <<
"gnorm";
952 os << std::setw(15) << std::left <<
"snorm";
953 os << std::setw(15) << std::left <<
"delta";
954 os << std::setw(10) << std::left <<
"#fval";
955 os << std::setw(10) << std::left <<
"#grad";
956 os << std::setw(10) << std::left <<
"#hess";
957 os << std::setw(10) << std::left <<
"#proj";
958 os << std::setw(10) << std::left <<
"tr_flag";
959 os << std::setw(10) << std::left <<
"iterSPG";
960 os << std::setw(10) << std::left <<
"flagSPG";
965 template<
typename Real>
967 std::ios_base::fmtflags osFlags(os.flags());
968 os << std::endl <<
"SPG Trust-Region Method (Type B, Bound Constraints)" << std::endl;
972 template<
typename Real>
974 std::ios_base::fmtflags osFlags(os.flags());
975 os << std::scientific << std::setprecision(6);
976 if ( state_->iter == 0 ) writeName(os);
977 if ( write_header ) writeHeader(os);
978 if ( state_->iter == 0 ) {
980 os << std::setw(6) << std::left << state_->iter;
981 os << std::setw(15) << std::left << state_->value;
982 os << std::setw(15) << std::left << state_->gnorm;
983 os << std::setw(15) << std::left <<
"---";
984 os << std::setw(15) << std::left << state_->searchSize;
985 os << std::setw(10) << std::left << state_->nfval;
986 os << std::setw(10) << std::left << state_->ngrad;
987 os << std::setw(10) << std::left << nhess_;
988 os << std::setw(10) << std::left << state_->nproj;
989 os << std::setw(10) << std::left <<
"---";
990 os << std::setw(10) << std::left <<
"---";
991 os << std::setw(10) << std::left <<
"---";
996 os << std::setw(6) << std::left << state_->iter;
997 os << std::setw(15) << std::left << state_->value;
998 os << std::setw(15) << std::left << state_->gnorm;
999 os << std::setw(15) << std::left << state_->snorm;
1000 os << std::setw(15) << std::left << state_->searchSize;
1001 os << std::setw(10) << std::left << state_->nfval;
1002 os << std::setw(10) << std::left << state_->ngrad;
1003 os << std::setw(10) << std::left << nhess_;
1004 os << std::setw(10) << std::left << state_->nproj;
1005 os << std::setw(10) << std::left << TRflag_;
1006 os << std::setw(10) << std::left << SPiter_;
1007 os << std::setw(10) << std::left << SPflag_;
Provides the interface to evaluate objective functions.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual void scale(const Real alpha)=0
Compute where .
void writeOutput(std::ostream &os, const bool write_header=false) const override
Print iterate status.
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual Real apply(const Vector< Real > &x) const
Apply to a dual vector. This is equivalent to the call .
virtual void plus(const Vector &x)=0
Compute , where .
void dpsg_simple(Vector< Real > &y, Real &q, Vector< Real > &gmod, const Vector< Real > &x, Real del, TrustRegionModel_U< Real > &model, Vector< Real > &pwa, Vector< Real > &pwa1, Vector< Real > &dwa, std::ostream &outStream=std::cout)
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
TrustRegionSPGAlgorithm(ParameterList &list, const Ptr< Secant< Real >> &secant=nullPtr)
void computeGradient(const Vector< Real > &x, Vector< Real > &g, Vector< Real > &pwa, Real del, Objective< Real > &obj, bool accept, Real >ol, Real &gnorm, std::ostream &outStream=std::cout) const
virtual void writeExitStatus(std::ostream &os) const
ESecant StringToESecant(std::string s)
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)
virtual void zero()
Set to zero vector.
Defines the linear algebra or vector space interface.
void dpsg(Vector< Real > &y, Real &q, Vector< Real > &gmod, const Vector< Real > &x, Real del, TrustRegionModel_U< Real > &model, Vector< Real > &ymin, Vector< Real > &pwa, Vector< Real > &pwa1, Vector< Real > &pwa2, Vector< Real > &pwa3, Vector< Real > &pwa4, Vector< Real > &pwa5, Vector< Real > &dwa, std::ostream &outStream=std::cout)
void writeName(std::ostream &os) const override
Print step name.
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
Real dgpstep(Vector< Real > &s, const Vector< Real > &w, const Vector< Real > &x, const Real alpha, std::ostream &outStream=std::cout) const
virtual Real dot(const Vector &x) const =0
Compute where .
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &s, Real &tol) override
Apply Hessian approximation to vector.
ETRFlag
Enumation of flags used by trust-region solvers.
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
std::string NumberToString(T Number)
Provides the interface to evaluate trust-region model functions.
void initialize(Vector< Real > &x, const Vector< Real > &g, Real ftol, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout)
void writeHeader(std::ostream &os) const override
Print iterate header.
Provides interface for and implements limited-memory secant operators.
void dproj(Vector< Real > &x, const Vector< Real > &x0, Real del, Vector< Real > &y0, Vector< Real > &y1, Vector< Real > &yc, Vector< Real > &pwa, std::ostream &outStream=std::cout) const
Provides an interface to check status of optimization algorithms.
std::string ETRFlagToString(ETRFlag trf)
Provides the interface to apply upper and lower bound constraints.
Real optimalityCriterion(const Vector< Real > &x, const Vector< Real > &g, Vector< Real > &primal, std::ostream &outStream=std::cout) const
void initialize(const Vector< Real > &x, const Vector< Real > &g)
virtual void set(const Vector &x)
Set where .
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...
virtual Real norm() const =0
Returns where .
Real computeValue(Real inTol, Real &outTol, Real pRed, Real &fold, int iter, const Vector< Real > &x, const Vector< Real > &xold, Objective< Real > &obj)