10#ifndef ROL_TYPEB_PRIMALDUALACTIVESETALGORITHM_DEF_HPP
11#define ROL_TYPEB_PRIMALDUALACTIVESETALGORITHM_DEF_HPP
16template<
typename Real>
26 secantName_ = list.sublist(
"General").sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS");
31 secantName_ = list.sublist(
"General").sublist(
"Secant").get(
"User Defined Secant Name",
32 "Unspecified User Defined Secant Method");
34 useSecantHessVec_ = list.sublist(
"General").sublist(
"Secant").get(
"Use as Hessian",
false);
35 useSecantPrecond_ = list.sublist(
"General").sublist(
"Secant").get(
"Use as Preconditioner",
false);
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);
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);
48 verbosity_ = list.sublist(
"General").get(
"Output Level", 0);
52template<
typename Real>
57 std::ostream &outStream) {
59 if (
proj_ == nullPtr) {
60 proj_ = makePtr<PolyhedralProjection<Real>>(makePtrFromRef(bnd));
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_);
70 krylov_ = makePtr<GMRES<Real>>(list);
82 proj_->project(x,outStream);
83 state_->iterateVec->set(x);
90 state_->stepVec->axpy(-one,x);
95template<
typename Real>
100 std::ostream &outStream ) {
101 const Real
zero(0), one(1);
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;
109 mu =
proj_->getMultiplier()->clone(); mu->zero();
110 b =
proj_->getResidual()->clone();
111 r =
proj_->getResidual()->clone();
113 lambda->set(
state_->gradientVec->dual());
117 Ptr<LinearOperator<Real>> hessian, precond;
127 proj_->getLinearConstraint()->value(*r,x,tol);
134 xlam->axpy(
scale_,*lambda);
141 xtmp->axpy(-one,*
state_->iterateVec);
146 xtmp->axpy(-one,*
state_->iterateVec);
159 proj_->getLinearConstraint()->applyJacobian(*b,*As,*
state_->iterateVec,tol);
163 gtmp->plus(*
state_->gradientVec);
170 rnorm = std::sqrt(gtmp->dot(*gtmp)+b->dot(*b));
173 rnorm = gtmp->norm();
175 if ( rnorm >
zero ) {
178 hessian = makePtr<HessianPDAS_Poly>(makePtrFromRef(obj),makePtrFromRef(bnd),
181 precond = makePtr<PrecondPDAS_Poly>(makePtrFromRef(obj),makePtrFromRef(bnd),
189 hessian = makePtr<HessianPDAS>(makePtrFromRef(obj),makePtrFromRef(bnd),
191 precond = makePtr<PrecondPDAS>(makePtrFromRef(obj),makePtrFromRef(bnd),
198 state_->stepVec->plus(*As);
204 snorm =
state_->stepVec->norm();
213 proj_->getLinearConstraint()->applyAdjointJacobian(*dwa,*mu,*
state_->iterateVec,tol);
216 gtmp->plus(*
state_->gradientVec);
218 lambda->set(gtmp->dual());
231 if ( snorm <
stol_*xnorm ) {
245 state_->iterateVec->set(x);
252 gtmp->set(*
state_->gradientVec);
255 xtmp->set(x); xtmp->axpy(-one,
state_->gradientVec->dual());
256 proj_->project(*xtmp,outStream);
258 state_->gnorm = xtmp->norm();
269template<
typename Real>
271 std::ios_base::fmtflags osFlags(os.flags());
273 os << std::string(114,
'-') << std::endl;
275 os <<
"Primal Dual Active Set Newton's Method";
278 os <<
"Primal Dual Active Set Quasi-Newton Method with " <<
secantName_ <<
" Hessian approximation";
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;
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;
293 os <<
" iterK - Number of Krylov iterations" << std::endl << std::endl;
294 os <<
" flagK - Krylov flag" << std::endl;
300 os <<
" feasible - Is iterate feasible?" << std::endl;
301 os << std::string(114,
'-') << std::endl;
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";
312 os << std::setw(10) << std::left <<
"iterPDAS";
313 os << std::setw(10) << std::left <<
"flagPDAS";
314 os << std::setw(10) << std::left <<
"iterK";
317 os << std::setw(10) << std::left <<
"iterK";
318 os << std::setw(10) << std::left <<
"flagK";
320 os << std::setw(10) << std::left <<
"feasible";
325template<
typename Real>
327 std::ios_base::fmtflags osFlags(os.flags());
329 os << std::endl <<
"Primal Dual Active Set Newton's Method (Type B, Bound Constraints)" << std::endl;
332 os << std::endl <<
"Primal Dual Active Set Quasi-Newton Method with "
333 <<
secantName_ <<
" Hessian approximation" << std::endl;
338template<
typename Real>
340 std::ios_base::fmtflags osFlags(os.flags());
341 os << std::scientific << std::setprecision(6);
344 if (
state_->iter == 0 ) {
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 <<
"---";
355 os << std::setw(10) << std::left <<
"---";
358 os << std::setw(10) << std::left <<
"YES";
361 os << std::setw(10) << std::left <<
"NO";
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;
374 os << std::setw(10) << std::left <<
iter_;
375 os << std::setw(10) << std::left <<
flag_;
383 os << std::setw(10) << std::left <<
"YES";
386 os << std::setw(10) << std::left <<
"NO";
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).
std::string krylovName_
Krylov name.
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.
int iter_
PDAS iteration counter.
ESecant esec_
Secant type.
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)
int flag_
PDAS termination flag.
Real neps_
-active set parameter
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...
std::string secantName_
Secant name.
int maxit_
Maximum number of PDAS steps (default: 10).
Ptr< Krylov< Real > > krylov_
Krylov solver object (used for inexact Newton).
Real itol_
Inexact Krylov tolerance.
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)
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
ESecant StringToESecant(std::string s)
ROL::Ptr< Secant< Real > > SecantFactory(ROL::ParameterList &parlist, ESecantMode mode=SECANTMODE_BOTH)
std::string ECGFlagToString(ECGFlag cgf)