ROL
ROL_TypeB_ColemanLiAlgorithm_Def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_TYPEB_COLEMANLIALGORITHM_DEF_HPP
11#define ROL_TYPEB_COLEMANLIALGORITHM_DEF_HPP
12
13namespace ROL {
14namespace TypeB {
15
16template<typename Real>
18 const Ptr<Secant<Real>> &secant) {
19 // Set status test
20 status_->reset();
21 status_->add(makePtr<StatusTest<Real>>(list));
22
23 ParameterList &trlist = list.sublist("Step").sublist("Trust Region");
24 // Trust-Region Parameters
25 state_->searchSize = trlist.get("Initial Radius", -1.0);
26 delMax_ = trlist.get("Maximum Radius", ROL_INF<Real>());
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);
36 // Nonmonotone Parameters
37 storageNM_ = trlist.get("Nonmonotone Storage Size", 0);
38 useNM_ = (storageNM_ <= 0 ? false : true);
39 // Krylov Parameters
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);
43 // Algorithm-Specific Parameters
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);
49 alphaMax_ = (alphaMax_ >= static_cast<Real>(1) ? static_cast<Real>(0.5) : alphaMax_);
50 // Output Parameters
51 verbosity_ = list.sublist("General").get("Output Level",0);
53 // Secant Information
54 useSecantPrecond_ = list.sublist("General").sublist("Secant").get("Use as Preconditioner", false);
55 useSecantHessVec_ = list.sublist("General").sublist("Secant").get("Use as Hessian", false);
59 // Initialize trust region model
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");
63 esec_ = StringToESecant(secantType);
64 }
65}
66
67template<typename Real>
69 const Vector<Real> &g,
70 Objective<Real> &obj,
72 std::ostream &outStream) {
73 const Real one(1);
74 hasEcon_ = true;
75 if (proj_ == nullPtr) {
76 proj_ = makePtr<PolyhedralProjection<Real>>(makePtrFromRef(bnd));
77 hasEcon_ = false;
78 }
79 // Initialize data
81 nhess_ = 0;
82 // Update approximate gradient and approximate objective function.
83 Real ftol = static_cast<Real>(0.1)*ROL_OVERFLOW<Real>();
84 proj_->getBoundConstraint()->projectInterior(x); state_->nproj++;
85 state_->iterateVec->set(x);
87 state_->value = obj.value(x,ftol);
88 state_->nfval++;
89 obj.gradient(*state_->gradientVec,x,ftol);
90 state_->ngrad++;
91 state_->stepVec->set(x);
92 state_->stepVec->axpy(-one,state_->gradientVec->dual());
93 proj_->project(*state_->stepVec,outStream); state_->nproj++;
94 state_->stepVec->axpy(-one,x);
95 state_->gnorm = state_->stepVec->norm();
96 state_->snorm = ROL_INF<Real>();
97 // Compute initial trust region radius if desired.
98 if ( state_->searchSize <= static_cast<Real>(0) ) {
99 state_->searchSize = state_->gradientVec->norm();
100 }
101 // Initialize null space projection
102 if (hasEcon_) {
103 rcon_ = makePtr<ReducedLinearConstraint<Real>>(proj_->getLinearConstraint(),
104 makePtrFromRef(bnd),
105 makePtrFromRef(x));
106 ns_ = makePtr<NullSpaceOperator<Real>>(rcon_,x,
107 *proj_->getResidual());
108 }
109}
110
111template<typename Real>
113 const Vector<Real> &g,
114 Objective<Real> &obj,
116 std::ostream &outStream ) {
117 const Real zero(0), one(1), half(0.5);
118 Real tol0 = std::sqrt(ROL_EPSILON<Real>());
119 Real tol(0), stol(0), snorm(0);
120 Real ftrial(0), pRed(0), rho(1), alpha(1);
121 // Initialize trust-region data
122 initialize(x,g,obj,bnd,outStream);
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();
126 // Initialize nonmonotone data
127 Real rhoNM(0), sigmac(0), sigmar(0), sBs(0), gs(0);
128 Real fr(state_->value), fc(state_->value), fmin(state_->value);
129 TRUtils::ETRFlag TRflagNM;
130 int L(0);
131
132 // Output
133 if (verbosity_ > 0) writeOutput(outStream,true);
134
135 while (status_->check(*state_)) {
136 // Build trust-region model (use only to encapsulate Hessian/secant)
137 model_->setData(obj,*state_->iterateVec,*state_->gradientVec,tol0);
138
139 // Run Truncated CG
140 // TODO: Model is: 1/2 (x-xk)' (B + Einv(D)) + g'(x-xk)
141 // applyHessian returns (B+Einv(D))v
142 SPflag_ = 0;
143 SPiter_ = 0;
144 tol = std::min(tol1_,tol2_*std::pow(state_->gnorm,spexp_));
145 stol = tol; //zero;
146 pwa5->set(state_->gradientVec->dual());
147 snorm = dtrpcg(*state_->stepVec,SPflag_,SPiter_,*state_->gradientVec,x,*pwa5,
148 state_->searchSize,*model_,bnd,tol,stol,
149 *pwa1,*dwa1,*pwa2,*dwa2,*pwa3,*pwa4,*dwa3,outStream);
150 if (verbosity_ > 1) {
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;
156 }
157
158 // Relax CG step so that it is interior
159 snorm = dgpstep(*pwa1,*state_->stepVec,x,one,outStream);
160 alpha = std::max(alphaMax_, one-snorm);
161 pwa1->scale(alpha);
162 state_->stepVec->set(*pwa1);
163 state_->snorm = alpha * snorm;
164 x.plus(*state_->stepVec);
165
166 // Compute predicted reduction
167 model_->hessVec(*dwa1,*pwa1,x,tol); nhess_++;
168 gs = state_->gradientVec->apply(*state_->stepVec);
169 sBs = dwa1->apply(*state_->stepVec);
170 pRed = - half * sBs - gs;
171
172 // Compute trial objective value
174 ftrial = obj.value(x,tol0);
175 state_->nfval++;
176
177 // Compute ratio of acutal and predicted reduction
179 TRUtils::analyzeRatio<Real>(rho,TRflag_,state_->value,ftrial,pRed,eps_,outStream,verbosity_>1);
180 if (useNM_) {
181 TRUtils::analyzeRatio<Real>(rhoNM,TRflagNM,fr,ftrial,pRed+sigmar,eps_,outStream,verbosity_>1);
182 TRflag_ = (rho < rhoNM ? TRflagNM : TRflag_);
183 rho = (rho < rhoNM ? rhoNM : rho );
184 }
185
186 // Update algorithm state
187 state_->iter++;
188 // Accept/reject step and update trust region radius
189 if ((rho < eta0_ && TRflag_ == TRUtils::SUCCESS) || (TRflag_ >= 2)) { // Step Rejected
190 x.set(*state_->iterateVec);
191 obj.update(x,UpdateType::Revert,state_->iter);
192 if (interpRad_ && (rho < zero && TRflag_ != TRUtils::TRNAN)) {
193 // Negative reduction, interpolate to find new trust-region radius
194 state_->searchSize = TRUtils::interpolateRadius<Real>(*state_->gradientVec,*state_->stepVec,
195 state_->snorm,pRed,state_->value,ftrial,state_->searchSize,gamma0_,gamma1_,eta2_,
196 outStream,verbosity_>1);
197 }
198 else { // Shrink trust-region radius
199 state_->searchSize = gamma1_*std::min(state_->snorm,state_->searchSize);
200 }
201 }
202 else if ((rho >= eta0_ && TRflag_ != TRUtils::NPOSPREDNEG)
203 || (TRflag_ == TRUtils::POSPREDNEG)) { // Step Accepted
204 state_->value = ftrial;
205 obj.update(x,UpdateType::Accept,state_->iter);
206 if (useNM_) {
207 sigmac += pRed; sigmar += pRed;
208 if (ftrial < fmin) { fmin = ftrial; fc = fmin; sigmac = zero; L = 0; }
209 else {
210 L++;
211 if (ftrial > fc) { fc = ftrial; sigmac = zero; }
212 if (L == storageNM_) { fr = fc; sigmar = sigmac; }
213 }
214 }
215 // Increase trust-region radius
216 if (rho >= eta2_) state_->searchSize = std::min(gamma2_*state_->searchSize, delMax_);
217 // Compute gradient at new iterate
218 dwa1->set(*state_->gradientVec);
219 obj.gradient(*state_->gradientVec,x,tol0);
220 state_->ngrad++;
221 state_->gnorm = TypeB::Algorithm<Real>::optimalityCriterion(x,*state_->gradientVec,*pwa1,outStream);
222 state_->iterateVec->set(x);
223 // Update secant information in trust-region model
224 model_->update(x,*state_->stepVec,*dwa1,*state_->gradientVec,
225 state_->snorm,state_->iter);
226 }
227
228 // Update Output
229 if (verbosity_ > 0) writeOutput(outStream,writeHeader_);
230 }
232}
233
234template<typename Real>
236 const Vector<Real> &x, const Real alpha,
237 std::ostream &outStream) const {
238 s.set(x); s.axpy(alpha,w);
239 proj_->getBoundConstraint()->projectInterior(s); state_->nproj++;
240 s.axpy(static_cast<Real>(-1),x);
241 return s.norm();
242}
243
244template<typename Real>
246 const Real ptp,
247 const Real ptx,
248 const Real del) const {
249 const Real zero(0);
250 Real dsq = del*del;
251 Real rad = ptx*ptx + ptp*(dsq-xtx);
252 rad = std::sqrt(std::max(rad,zero));
253 Real sigma(0);
254 if (ptx > zero) {
255 sigma = (dsq-xtx)/(ptx+rad);
256 }
257 else if (rad > zero) {
258 sigma = (rad-ptx)/ptp;
259 }
260 else {
261 sigma = zero;
262 }
263 return sigma;
264}
265
266template<typename Real>
267Real ColemanLiAlgorithm<Real>::dtrpcg(Vector<Real> &w, int &iflag, int &iter,
268 const Vector<Real> &g, const Vector<Real> &x,
269 const Vector<Real> &gdual,
270 const Real del, TrustRegionModel_U<Real> &model,
272 const Real tol, const Real stol,
274 Vector<Real> &t, Vector<Real> &pwa1,
275 Vector<Real> &pwa2, Vector<Real> &dwa,
276 std::ostream &outStream) const {
277 // p = step (primal)
278 // q = hessian applied to step p (dual)
279 // t = gradient (dual)
280 // r = preconditioned gradient (primal)
281 Real tol0 = std::sqrt(ROL_EPSILON<Real>());
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);
285 iter = 0; iflag = 0;
286 // Initialize step
287 w.zero();
288 // Compute residual
289 t.set(g); t.scale(-one);
290 // Preconditioned residual
291 applyPrecond(r,t,x,gdual,model,bnd,tol0,dwa,pwa1);
292 //rho = r.dot(t.dual());
293 rho = r.apply(t);
294 // Initialize direction
295 p.set(r);
296 pMp = (!hasEcon_ ? rho : p.dot(p)); // If no equality constraint, used preconditioned norm
297 // Iterate CG
298 for (iter = 0; iter < maxit_; ++iter) {
299 // Apply Hessian to direction dir
300 applyHessian(q,p,x,gdual,model,bnd,tol0,pwa1,pwa2);
301 // Compute sigma such that ||s+sigma*dir|| = del
302 //kappa = p.dot(q.dual());
303 kappa = p.apply(q);
304 alpha = (kappa>zero) ? rho/kappa : zero;
305 sigma = dtrqsol(sMs,pMp,sMp,del);
306 // Check for negative curvature or if iterate exceeds trust region
307 if (kappa <= zero || alpha >= sigma) {
308 w.axpy(sigma,p);
309 sMs = del*del;
310 iflag = (kappa<=zero) ? 2 : 3;
311 break;
312 }
313 // Update iterate and residuals
314 w.axpy(alpha,p);
315 t.axpy(-alpha,q);
316 applyPrecond(r,t,x,gdual,model,bnd,tol0,dwa,pwa1);
317 // Exit if residual tolerance is met
318 //rtr = r.dot(t.dual());
319 rtr = r.apply(t);
320 tnorm = t.norm();
321 if (rtr <= stol*stol || tnorm <= tol) {
322 sMs = sMs + two*alpha*sMp + alpha*alpha*pMp;
323 iflag = 0;
324 break;
325 }
326 // Compute p = r + beta * p
327 beta = rtr/rho;
328 p.scale(beta); p.plus(r);
329 rho = rtr;
330 // Update dot products
331 // sMs = <s, inv(M)s> or <s, s> if equality constraint present
332 // sMp = <s, inv(M)p> or <s, p> if equality constraint present
333 // pMp = <p, inv(M)p> or <p, p> if equality constraint present
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;
337 }
338 // Check iteration count
339 if (iter == maxit_) {
340 iflag = 1;
341 }
342 if (iflag != 1) {
343 iter++;
344 }
345 return std::sqrt(sMs); // w.norm();
346}
347
348template<typename Real>
350 const Vector<Real> &v,
351 const Vector<Real> &x,
352 const Vector<Real> &g,
354 Vector<Real> &pwa) const {
355 bnd.applyInverseScalingFunction(pwa,v,x,g);
356 bnd.applyScalingFunctionJacobian(Cv,pwa,x,g);
357}
358
359template<typename Real>
361 const Vector<Real> &v,
362 const Vector<Real> &x,
363 const Vector<Real> &g,
366 Real &tol,
367 Vector<Real> &pwa1,
368 Vector<Real> &pwa2) const {
369 model.hessVec(hv,v,x,tol); nhess_++;
370 applyC(pwa2,v,x,g,bnd,pwa1);
371 hv.plus(pwa2.dual());
372}
373
374template<typename Real>
376 const Vector<Real> &v,
377 const Vector<Real> &x,
378 const Vector<Real> &g,
381 Real &tol,
382 Vector<Real> &dwa,
383 Vector<Real> &pwa) const {
384 model.precond(hv,v,x,tol);
385}
386
387template<typename Real>
388void ColemanLiAlgorithm<Real>::writeHeader( std::ostream& os ) const {
389 std::ios_base::fmtflags osFlags(os.flags());
390 if (verbosity_ > 1) {
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;
402 os << std::endl;
403 os << " tr_flag - Trust-Region flag" << std::endl;
404 for( int flag = TRUtils::SUCCESS; flag != TRUtils::UNDEFINED; ++flag ) {
405 os << " " << NumberToString(flag) << " - "
406 << TRUtils::ETRFlagToString(static_cast<TRUtils::ETRFlag>(flag)) << std::endl;
407 }
408 os << std::endl;
409 os << " iterCG - Number of Truncated CG iterations" << std::endl << std::endl;
410 os << " flagGC - Trust-Region Truncated CG flag" << std::endl;
411 for( int flag = CG_FLAG_SUCCESS; flag != CG_FLAG_UNDEFINED; ++flag ) {
412 os << " " << NumberToString(flag) << " - "
413 << ECGFlagToString(static_cast<ECGFlag>(flag)) << std::endl;
414 }
415 os << std::string(114,'-') << std::endl;
416 }
417 os << " ";
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";
430 os << std::endl;
431 os.flags(osFlags);
432}
433
434template<typename Real>
435void ColemanLiAlgorithm<Real>::writeName( std::ostream& os ) const {
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;
438 os.flags(osFlags);
439}
440
441template<typename Real>
442void ColemanLiAlgorithm<Real>::writeOutput( std::ostream& os, bool write_header ) const {
443 std::ios_base::fmtflags osFlags(os.flags());
444 os << std::scientific << std::setprecision(6);
445 if ( state_->iter == 0 ) writeName(os);
446 if ( write_header ) writeHeader(os);
447 if ( state_->iter == 0 ) {
448 os << " ";
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 << "---";
461 os << std::endl;
462 }
463 else {
464 os << " ";
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_;
477 os << std::endl;
478 }
479 os.flags(osFlags);
480}
481
482} // namespace TypeB
483} // namespace ROL
484
485#endif
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)
Definition ROL_Types.hpp:47
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
Definition ROL_Types.hpp:57
ESecant StringToESecant(std::string s)
ESecantMode
@ SECANTMODE_FORWARD
@ SECANTMODE_INVERSE
@ SECANTMODE_BOTH
@ CG_FLAG_UNDEFINED
@ CG_FLAG_SUCCESS
Real ROL_OVERFLOW(void)
Platform-dependent maximum double.
Definition ROL_Types.hpp:68
std::string ECGFlagToString(ECGFlag cgf)
Real ROL_INF(void)
Definition ROL_Types.hpp:71