ROL
ROL_TypeB_PrimalDualActiveSetAlgorithm_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_PRIMALDUALACTIVESETALGORITHM_DEF_HPP
11#define ROL_TYPEB_PRIMALDUALACTIVESETALGORITHM_DEF_HPP
12
13namespace ROL {
14namespace TypeB {
15
16template<typename Real>
18 const Ptr<Secant<Real>> &secant)
20 neps_(-ROL_EPSILON<Real>()), itol_(std::sqrt(ROL_EPSILON<Real>())), hasPoly_(true) {
21 // Set status test
22 status_->reset();
23 status_->add(makePtr<StatusTest<Real>>(list));
24
25 if ( secant_ == nullPtr ) {
26 secantName_ = list.sublist("General").sublist("Secant").get("Type","Limited-Memory BFGS");
29 }
30 else {
31 secantName_ = list.sublist("General").sublist("Secant").get("User Defined Secant Name",
32 "Unspecified User Defined Secant Method");
33 }
34 useSecantHessVec_ = list.sublist("General").sublist("Secant").get("Use as Hessian", false);
35 useSecantPrecond_ = list.sublist("General").sublist("Secant").get("Use as Preconditioner", false);
36
37 // Algorithmic parameters
38 maxit_ = list.sublist("Step").sublist("Primal Dual Active Set").get("Iteration Limit", 10);
39 stol_ = list.sublist("Step").sublist("Primal Dual Active Set").get("Relative Step Tolerance", 1e-8);
40 gtol_ = list.sublist("Step").sublist("Primal Dual Active Set").get("Relative Gradient Tolerance", 1e-6);
41 scale_ = list.sublist("Step").sublist("Primal Dual Active Set").get("Dual Scaling", 1.0);
42
43 // Krylov parameters
44 atolKrylov_ = list.sublist("General").sublist("Krylov").get("Absolute Tolerance", 1e-4);
45 rtolKrylov_ = list.sublist("General").sublist("Krylov").get("Relative Tolerance", 1e-2);
46 maxitKrylov_ = list.sublist("General").sublist("Krylov").get("Iteration Limit", 100);
47
48 verbosity_ = list.sublist("General").get("Output Level", 0);
50}
51
52template<typename Real>
54 const Vector<Real> &g,
55 Objective<Real> &obj,
57 std::ostream &outStream) {
58 // Initialize projection operator
59 if (proj_ == nullPtr) {
60 proj_ = makePtr<PolyhedralProjection<Real>>(makePtrFromRef(bnd));
61 hasPoly_ = false;
62 }
63 // Create Krylov solver
64 if (hasPoly_) {
65 ParameterList list;
66 list.sublist("General").sublist("Krylov").set("Absolute Tolerance", atolKrylov_);
67 list.sublist("General").sublist("Krylov").set("Relative Tolerance", rtolKrylov_);
68 list.sublist("General").sublist("Krylov").set("Iteration Limit", maxitKrylov_);
69 krylovName_ = "GMRES";
70 krylov_ = makePtr<GMRES<Real>>(list);
71 }
72 else {
73 krylovName_ = "CR";
74 krylov_ = makePtr<ConjugateResiduals<Real>>(atolKrylov_,rtolKrylov_,maxitKrylov_);
75 }
77 // Initialize data
78 const Real one(1);
80 // Update approximate gradient and approximate objective function.
81 Real ftol = std::sqrt(ROL_EPSILON<Real>());
82 proj_->project(x,outStream);
83 state_->iterateVec->set(x);
85 state_->value = obj.value(x,ftol); state_->nfval++;
86 obj.gradient(*state_->gradientVec,x,ftol); state_->ngrad++;
87 state_->stepVec->set(x);
88 state_->stepVec->axpy(-one,state_->gradientVec->dual());
89 proj_->project(*state_->stepVec,outStream);
90 state_->stepVec->axpy(-one,x);
91 state_->gnorm = state_->stepVec->norm();
92 state_->snorm = ROL_INF<Real>();
93}
94
95template<typename Real>
97 const Vector<Real> &g,
98 Objective<Real> &obj,
100 std::ostream &outStream ) {
101 const Real zero(0), one(1);
102 // Initialize PDAS data
103 initialize(x,g,obj,bnd,outStream);
104 Ptr<Vector<Real>> xlam = x.clone(), xtmp = x.clone(), As = x.clone();
105 Ptr<Vector<Real>> lambda = x.clone(), gtmp = g.clone();
106 Ptr<Vector<Real>> pwa = x.clone(), dwa = g.clone();
107 Ptr<Vector<Real>> mu, b, r;
108 if (hasPoly_) {
109 mu = proj_->getMultiplier()->clone(); mu->zero();
110 b = proj_->getResidual()->clone();
111 r = proj_->getResidual()->clone();
112 }
113 lambda->set(state_->gradientVec->dual());
114 lambda->scale(-one);
115 Real xnorm(0), snorm(0), rnorm(0), tol(std::sqrt(ROL_EPSILON<Real>()));
116
117 Ptr<LinearOperator<Real>> hessian, precond;
118
119 // Output
120 if (verbosity_ > 0) writeOutput(outStream,true);
121
122 while (status_->check(*state_)) {
123 totalKrylov_ = 0;
124 state_->stepVec->zero();
125 xnorm = x.norm();
126 if (hasPoly_) {
127 proj_->getLinearConstraint()->value(*r,x,tol);
128 }
129 for ( iter_ = 0; iter_ < maxit_; iter_++ ) {
130 /********************************************************************/
131 // MODIFY ITERATE VECTOR TO CHECK ACTIVE SET
132 /********************************************************************/
133 xlam->set(x); // xlam = x0
134 xlam->axpy(scale_,*lambda); // xlam = x0 + c*lambda
135 /********************************************************************/
136 // PROJECT x ONTO PRIMAL DUAL FEASIBLE SET
137 /********************************************************************/
138 As->zero(); // As = 0
139
140 xtmp->set(*bnd.getUpperBound()); // xtmp = u
141 xtmp->axpy(-one,*state_->iterateVec); // xtmp = u - x
142 bnd.pruneUpperInactive(*xtmp,*xlam,neps_); // xtmp = A(u - x)
143 As->plus(*xtmp); // As += A(u - x)
144
145 xtmp->set(*bnd.getLowerBound()); // xtmp = l
146 xtmp->axpy(-one,*state_->iterateVec); // xtmp = l - x
147 bnd.pruneLowerInactive(*xtmp,*xlam,neps_); // xtmp = A(l - x)
148 As->plus(*xtmp); // As += A(l - x)
149 /********************************************************************/
150 // APPLY HESSIAN TO ACTIVE COMPONENTS OF s AND REMOVE INACTIVE
151 /********************************************************************/
152 if ( useSecantHessVec_ && secant_ != nullPtr ) { // gtmp = H*As
153 secant_->applyB(*gtmp,*As);
154 }
155 else {
156 obj.hessVec(*gtmp,*As,*state_->iterateVec,itol_);
157 }
158 if (hasPoly_) {
159 proj_->getLinearConstraint()->applyJacobian(*b,*As,*state_->iterateVec,tol);
160 b->plus(*r);
161 b->scale(-one);
162 }
163 gtmp->plus(*state_->gradientVec); // gtmp = g + H*As + ...
164 gtmp->scale(-one); // gtmp = -(g + H*As + ...)
165 bnd.pruneActive(*gtmp,*xlam,neps_); // gtmp = -I(g + H*As + ...)
166 /********************************************************************/
167 // SOLVE REDUCED NEWTON SYSTEM
168 /********************************************************************/
169 if (hasPoly_) {
170 rnorm = std::sqrt(gtmp->dot(*gtmp)+b->dot(*b));
171 }
172 else {
173 rnorm = gtmp->norm();
174 }
175 if ( rnorm > zero ) {
176 if (hasPoly_) {
177 // Initialize Hessian and preconditioner
178 hessian = makePtr<HessianPDAS_Poly>(makePtrFromRef(obj),makePtrFromRef(bnd),
179 proj_->getLinearConstraint(),state_->iterateVec,xlam,neps_,secant_,
180 useSecantHessVec_,pwa,dwa);
181 precond = makePtr<PrecondPDAS_Poly>(makePtrFromRef(obj),makePtrFromRef(bnd),
182 state_->iterateVec,xlam,neps_,secant_,useSecantPrecond_,dwa);
183 PartitionedVector<Real> rhs(std::vector<Ptr<Vector<Real>>>({gtmp,b}));
184 PartitionedVector<Real> sol(std::vector<Ptr<Vector<Real>>>({state_->stepVec,mu}));
185 krylov_->run(sol,*hessian,rhs,*precond,iterKrylov_,flagKrylov_);
186 }
187 else {
188 // Initialize Hessian and preconditioner
189 hessian = makePtr<HessianPDAS>(makePtrFromRef(obj),makePtrFromRef(bnd),
190 state_->iterateVec,xlam,neps_,secant_,useSecantHessVec_,pwa);
191 precond = makePtr<PrecondPDAS>(makePtrFromRef(obj),makePtrFromRef(bnd),
192 state_->iterateVec,xlam,neps_,secant_,useSecantPrecond_,dwa);
193 krylov_->run(*state_->stepVec,*hessian,*gtmp,*precond,iterKrylov_,flagKrylov_);
194 }
196 bnd.pruneActive(*state_->stepVec,*xlam,neps_); // s <- Is
197 }
198 state_->stepVec->plus(*As); // s = Is + As
199 /********************************************************************/
200 // UPDATE STEP AND MULTIPLIER
201 /********************************************************************/
202 x.set(*state_->iterateVec);
203 x.plus(*state_->stepVec);
204 snorm = state_->stepVec->norm();
205 // Compute gradient of Lagrangian for QP
206 if ( useSecantHessVec_ && secant_ != nullPtr ) {
207 secant_->applyB(*gtmp,*state_->stepVec);
208 }
209 else {
210 obj.hessVec(*gtmp,*state_->stepVec,*state_->iterateVec,itol_);
211 }
212 if (hasPoly_) {
213 proj_->getLinearConstraint()->applyAdjointJacobian(*dwa,*mu,*state_->iterateVec,tol);
214 gtmp->plus(*dwa);
215 }
216 gtmp->plus(*state_->gradientVec);
217 gtmp->scale(-one);
218 lambda->set(gtmp->dual());
219 // Compute criticality measure
220 xtmp->set(x);
221 xtmp->plus(*lambda);
222 bnd.project(*xtmp);
223 xtmp->axpy(-one,x);
224 // Update multiplier
225 bnd.pruneInactive(*lambda,*xlam,neps_);
226 // Check stopping conditions
227 if ( xtmp->norm() < gtol_*state_->gnorm ) {
228 flag_ = 0;
229 break;
230 }
231 if ( snorm < stol_*xnorm ) {
232 flag_ = 2;
233 break;
234 }
235 }
236 if ( iter_ == maxit_ ) {
237 flag_ = 1;
238 }
239 else {
240 iter_++;
241 }
242
243 // Update algorithm state
244 state_->iter++;
245 state_->iterateVec->set(x);
246 feasible_ = bnd.isFeasible(x);
247 state_->snorm = snorm;
248 obj.update(x,UpdateType::Accept,state_->iter);
249 state_->value = obj.value(x,tol); state_->nfval++;
250
251 if ( secant_ != nullPtr ) {
252 gtmp->set(*state_->gradientVec);
253 }
254 obj.gradient(*state_->gradientVec,x,tol); state_->ngrad++;
255 xtmp->set(x); xtmp->axpy(-one,state_->gradientVec->dual());
256 proj_->project(*xtmp,outStream);
257 xtmp->axpy(-one,x);
258 state_->gnorm = xtmp->norm();
259 if ( secant_ != nullPtr ) {
260 secant_->updateStorage(x,*state_->gradientVec,*gtmp,*state_->stepVec,state_->snorm,state_->iter+1);
261 }
262
263 // Update Output
264 if (verbosity_ > 0) writeOutput(outStream,writeHeader_);
265 }
267}
268
269template<typename Real>
271 std::ios_base::fmtflags osFlags(os.flags());
272 if (verbosity_ > 1) {
273 os << std::string(114,'-') << std::endl;
274 if (!useSecantHessVec_) {
275 os << "Primal Dual Active Set Newton's Method";
276 }
277 else {
278 os << "Primal Dual Active Set Quasi-Newton Method with " << secantName_ << " Hessian approximation";
279 }
280 os << " status output definitions" << std::endl << std::endl;
281 os << " iter - Number of iterates (steps taken)" << std::endl;
282 os << " value - Objective function value" << std::endl;
283 os << " gnorm - Norm of the gradient" << std::endl;
284 os << " snorm - Norm of the step (update to optimization vector)" << std::endl;
285 os << " #fval - Cumulative number of times the objective function was evaluated" << std::endl;
286 os << " #grad - Cumulative number of times the gradient was computed" << std::endl;
287 if (maxit_ > 1) {
288 os << " iterPDAS - Number of Primal Dual Active Set iterations" << std::endl << std::endl;
289 os << " flagPDAS - Primal Dual Active Set flag" << std::endl;
290 os << " iterK - Number of Krylov iterations" << std::endl << std::endl;
291 }
292 else {
293 os << " iterK - Number of Krylov iterations" << std::endl << std::endl;
294 os << " flagK - Krylov flag" << std::endl;
295 for( int flag = CG_FLAG_SUCCESS; flag != CG_FLAG_UNDEFINED; ++flag ) {
296 os << " " << NumberToString(flag) << " - "
297 << ECGFlagToString(static_cast<ECGFlag>(flag)) << std::endl;
298 }
299 }
300 os << " feasible - Is iterate feasible?" << std::endl;
301 os << std::string(114,'-') << std::endl;
302 }
303
304 os << " ";
305 os << std::setw(6) << std::left << "iter";
306 os << std::setw(15) << std::left << "value";
307 os << std::setw(15) << std::left << "gnorm";
308 os << std::setw(15) << std::left << "snorm";
309 os << std::setw(10) << std::left << "#fval";
310 os << std::setw(10) << std::left << "#grad";
311 if (maxit_ > 1) {
312 os << std::setw(10) << std::left << "iterPDAS";
313 os << std::setw(10) << std::left << "flagPDAS";
314 os << std::setw(10) << std::left << "iterK";
315 }
316 else {
317 os << std::setw(10) << std::left << "iterK";
318 os << std::setw(10) << std::left << "flagK";
319 }
320 os << std::setw(10) << std::left << "feasible";
321 os << std::endl;
322 os.flags(osFlags);
323}
324
325template<typename Real>
326void PrimalDualActiveSetAlgorithm<Real>::writeName( std::ostream& os ) const {
327 std::ios_base::fmtflags osFlags(os.flags());
328 if (!useSecantHessVec_) {
329 os << std::endl << "Primal Dual Active Set Newton's Method (Type B, Bound Constraints)" << std::endl;
330 }
331 else {
332 os << std::endl << "Primal Dual Active Set Quasi-Newton Method with "
333 << secantName_ << " Hessian approximation" << std::endl;
334 }
335 os.flags(osFlags);
336}
337
338template<typename Real>
339void PrimalDualActiveSetAlgorithm<Real>::writeOutput( std::ostream& os, bool write_header ) const {
340 std::ios_base::fmtflags osFlags(os.flags());
341 os << std::scientific << std::setprecision(6);
342 if ( state_->iter == 0 ) writeName(os);
343 if ( write_header ) writeHeader(os);
344 if ( state_->iter == 0 ) {
345 os << " ";
346 os << std::setw(6) << std::left << state_->iter;
347 os << std::setw(15) << std::left << state_->value;
348 os << std::setw(15) << std::left << state_->gnorm;
349 os << std::setw(15) << std::left << "---";
350 os << std::setw(10) << std::left << state_->nfval;
351 os << std::setw(10) << std::left << state_->ngrad;
352 os << std::setw(10) << std::left << "---";
353 os << std::setw(10) << std::left << "---";
354 if (maxit_ > 1) {
355 os << std::setw(10) << std::left << "---";
356 }
357 if ( feasible_ ) {
358 os << std::setw(10) << std::left << "YES";
359 }
360 else {
361 os << std::setw(10) << std::left << "NO";
362 }
363 os << std::endl;
364 }
365 else {
366 os << " ";
367 os << std::setw(6) << std::left << state_->iter;
368 os << std::setw(15) << std::left << state_->value;
369 os << std::setw(15) << std::left << state_->gnorm;
370 os << std::setw(15) << std::left << state_->snorm;
371 os << std::setw(10) << std::left << state_->nfval;
372 os << std::setw(10) << std::left << state_->ngrad;
373 if (maxit_ > 1) {
374 os << std::setw(10) << std::left << iter_;
375 os << std::setw(10) << std::left << flag_;
376 os << std::setw(10) << std::left << totalKrylov_;
377 }
378 else {
379 os << std::setw(10) << std::left << iterKrylov_;
380 os << std::setw(10) << std::left << flagKrylov_;
381 }
382 if ( feasible_ ) {
383 os << std::setw(10) << std::left << "YES";
384 }
385 else {
386 os << std::setw(10) << std::left << "NO";
387 }
388 os << std::endl;
389 }
390 os.flags(osFlags);
391}
392
393} // namespace TypeB
394} // namespace ROL
395
396#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 const Ptr< const Vector< Real > > getLowerBound(void) const
Return the ref count pointer to the lower bound vector.
virtual bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
void pruneLowerInactive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the -inactive set.
void pruneInactive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the -inactive set.
void pruneUpperInactive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the -inactive set.
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=Real(0))
Set variables to zero if they correspond to the -active set.
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
virtual const Ptr< const Vector< Real > > getUpperBound(void) const
Return the ref count pointer to the upper bound vector.
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 hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
virtual void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
Defines the linear algebra of vector space on a generic partitioned vector.
Provides interface for and implements limited-memory secant operators.
Provides an interface to check status of optimization algorithms.
Ptr< PolyhedralProjection< Real > > proj_
void initialize(const Vector< Real > &x, const Vector< Real > &g)
virtual void writeExitStatus(std::ostream &os) const
const Ptr< AlgorithmState< Real > > state_
const Ptr< CombinedStatusTest< Real > > status_
int flagKrylov_
Termination flag for Krylov method (used for inexact Newton).
int iterKrylov_
Number of Krylov iterations (used for inexact Newton).
Real stol_
PDAS minimum step size stopping tolerance (default: 1e-8).
void writeName(std::ostream &os) const override
Print step name.
Real scale_
Scale for dual variables in the active set, (default: 1).
bool feasible_
Flag whether the current iterate is feasible or not.
PrimalDualActiveSetAlgorithm(ParameterList &list, const Ptr< Secant< Real > > &secant=nullPtr)
Ptr< Secant< Real > > secant_
Secant object (used for quasi-Newton).
Real atolKrylov_
Absolute tolerance for Krylov solve (default: 1e-4).
int totalKrylov_
Total number of Krylov iterations per PDAS iteration.
Real gtol_
PDAS gradient stopping tolerance (default: 1e-6).
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, std::ostream &outStream=std::cout)
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 maxit_
Maximum number of PDAS steps (default: 10).
Ptr< Krylov< Real > > krylov_
Krylov solver object (used for inexact Newton).
void writeHeader(std::ostream &os) const override
Print iterate header.
int maxitKrylov_
Maximum number of Krylov iterations (default: 100).
bool useSecantHessVec_
Whether or not to use to a secant approximation as the Hessian.
void writeOutput(std::ostream &os, const bool write_header=false) const override
Print iterate status.
Real rtolKrylov_
Relative tolerance for Krylov solve (default: 1e-2).
bool useSecantPrecond_
Whether or not to use a secant approximation to precondition inexact Newton.
Defines the linear algebra or vector space interface.
virtual Real norm() const =0
Returns where .
virtual void set(const Vector &x)
Set where .
virtual void plus(const Vector &x)=0
Compute , where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
EKrylov StringToEKrylov(std::string s)
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)
@ SECANT_USERDEFINED
@ CG_FLAG_UNDEFINED
@ CG_FLAG_SUCCESS
ROL::Ptr< Secant< Real > > SecantFactory(ROL::ParameterList &parlist, ESecantMode mode=SECANTMODE_BOTH)
std::string ECGFlagToString(ECGFlag cgf)
Real ROL_INF(void)
Definition ROL_Types.hpp:71