44 #ifndef ROL_PRIMALDUALACTIVESETSTEP_H
45 #define ROL_PRIMALDUALACTIVESETSTEP_H
54 #include "ROL_ParameterList.hpp"
132 template <
class Real>
156 ROL::Ptr<Vector<Real> >
x0_;
158 ROL::Ptr<Vector<Real> >
As_;
161 ROL::Ptr<Vector<Real> >
Ag_;
173 const ROL::Ptr<Objective<Real> >
obj_;
174 const ROL::Ptr<BoundConstraint<Real> >
bnd_;
175 const ROL::Ptr<Vector<Real> >
x_;
176 const ROL::Ptr<Vector<Real> >
xlam_;
177 ROL::Ptr<Vector<Real> >
v_;
188 const bool useSecant =
false )
192 if ( !useSecant || secant == ROL::nullPtr ) {
211 const ROL::Ptr<Objective<Real> >
obj_;
212 const ROL::Ptr<BoundConstraint<Real> >
bnd_;
213 const ROL::Ptr<Vector<Real> >
x_;
214 const ROL::Ptr<Vector<Real> >
xlam_;
215 ROL::Ptr<Vector<Real> >
v_;
226 const bool useSecant =
false )
229 v_ =
x_->dual().clone();
230 if ( !useSecant || secant == ROL::nullPtr ) {
265 obj.
gradient(*(step_state->gradientVec),x,tol);
267 xtmp_->axpy(-one,(step_state->gradientVec)->dual());
270 return xtmp_->norm();
287 res_(ROL::nullPtr),
Ag_(ROL::nullPtr),
rtmp_(ROL::nullPtr),
291 Real one(1), oem6(1.e-6), oem8(1.e-8);
293 maxit_ = parlist.sublist(
"Step").sublist(
"Primal Dual Active Set").get(
"Iteration Limit",10);
294 stol_ = parlist.sublist(
"Step").sublist(
"Primal Dual Active Set").get(
"Relative Step Tolerance",oem8);
295 gtol_ = parlist.sublist(
"Step").sublist(
"Primal Dual Active Set").get(
"Relative Gradient Tolerance",oem6);
296 scale_ = parlist.sublist(
"Step").sublist(
"Primal Dual Active Set").get(
"Dual Scaling", one);
298 esec_ =
StringToESecant(parlist.sublist(
"General").sublist(
"Secant").get(
"Type",
"Limited-Memory BFGS"));
299 useSecantHessVec_ = parlist.sublist(
"General").sublist(
"Secant").get(
"Use as Hessian",
false);
300 useSecantPrecond_ = parlist.sublist(
"General").sublist(
"Secant").get(
"Use as Preconditioner",
false);
302 secant_ = SecantFactory<Real>(parlist);
305 krylov_ = KrylovFactory<Real>(parlist);
324 Real
zero(0), one(1);
326 step_state->descentVec = s.
clone();
327 step_state->gradientVec = g.
clone();
328 step_state->searchSize =
zero;
342 Real tol = std::sqrt(ROL_EPSILON<Real>());
350 lambda_->set((step_state->gradientVec)->dual());
383 Real
zero(0), one(1);
386 res_->set(*(step_state->gradientVec));
414 itol_ = std::sqrt(ROL_EPSILON<Real>());
425 rtmp_->set(*(step_state->gradientVec));
428 Ag_->set(*(step_state->gradientVec));
438 ROL::Ptr<Objective<Real> > obj_ptr = ROL::makePtrFromRef(obj);
439 ROL::Ptr<BoundConstraint<Real> > con_ptr = ROL::makePtrFromRef(con);
440 ROL::Ptr<LinearOperator<Real> > hessian
441 = ROL::makePtr<HessianPD>(obj_ptr,con_ptr,
443 ROL::Ptr<LinearOperator<Real> > precond
444 = ROL::makePtr<PrecondPD>(obj_ptr,con_ptr,
471 res_->set(*(step_state->gradientVec));
493 if (
iter_ == maxit_ ) {
523 Real tol = std::sqrt(ROL_EPSILON<Real>());
528 if (
secant_ != ROL::nullPtr ) {
529 gtmp_->set(*(step_state->gradientVec));
534 if (
secant_ != ROL::nullPtr ) {
546 std::stringstream hist;
548 hist << std::setw(6) << std::left <<
"iter";
549 hist << std::setw(15) << std::left <<
"value";
550 hist << std::setw(15) << std::left <<
"gnorm";
551 hist << std::setw(15) << std::left <<
"snorm";
552 hist << std::setw(10) << std::left <<
"#fval";
553 hist << std::setw(10) << std::left <<
"#grad";
555 hist << std::setw(10) << std::left <<
"iterPDAS";
556 hist << std::setw(10) << std::left <<
"flagPDAS";
559 hist << std::setw(10) << std::left <<
"iterCR";
560 hist << std::setw(10) << std::left <<
"flagCR";
562 hist << std::setw(10) << std::left <<
"feasible";
573 std::stringstream hist;
574 hist <<
"\nPrimal Dual Active Set Newton's Method\n";
586 std::stringstream hist;
587 hist << std::scientific << std::setprecision(6);
588 if ( algo_state.
iter == 0 ) {
591 if ( print_header ) {
594 if ( algo_state.
iter == 0 ) {
596 hist << std::setw(6) << std::left << algo_state.
iter;
597 hist << std::setw(15) << std::left << algo_state.
value;
598 hist << std::setw(15) << std::left << algo_state.
gnorm;
603 hist << std::setw(6) << std::left << algo_state.
iter;
604 hist << std::setw(15) << std::left << algo_state.
value;
605 hist << std::setw(15) << std::left << algo_state.
gnorm;
606 hist << std::setw(15) << std::left << algo_state.
snorm;
607 hist << std::setw(10) << std::left << algo_state.
nfval;
608 hist << std::setw(10) << std::left << algo_state.
ngrad;
610 hist << std::setw(10) << std::left <<
iter_;
611 hist << std::setw(10) << std::left <<
flag_;
614 hist << std::setw(10) << std::left <<
iterCR_;
615 hist << std::setw(10) << std::left <<
flagCR_;
618 hist << std::setw(10) << std::left <<
"YES";
621 hist << std::setw(10) << std::left <<
"NO";
ROL::Ptr< Vector< Real > > v_
Implements the computation of optimization steps with the Newton primal-dual active set method...
Provides the interface to evaluate objective functions.
virtual const ROL::Ptr< const Vector< Real > > getUpperBound(void) const
Return the ref count pointer to the upper bound vector.
std::string printName(void) const
Print step name.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void plus(const Vector &x)=0
Compute , where .
std::string printHeader(void) const
Print iterate header.
int iter_
PDAS iteration counter.
int iterCR_
CR iteration counter.
virtual Real value(const Vector< Real > &x, Real &tol)=0
Compute value.
Provides the interface to compute optimization steps.
const ROL::Ptr< BoundConstraint< Real > > bnd_
virtual void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
const ROL::Ptr< Vector< Real > > xlam_
Contains definitions of custom data types in ROL.
Real neps_
-active set parameter
const ROL::Ptr< Secant< Real > > secant_
HessianPD(const ROL::Ptr< Objective< Real > > &obj, const ROL::Ptr< BoundConstraint< Real > > &bnd, const ROL::Ptr< Vector< Real > > &x, const ROL::Ptr< Vector< Real > > &xlam, const Real eps=0, const ROL::Ptr< Secant< Real > > &secant=ROL::nullPtr, const bool useSecant=false)
void pruneActive(Vector< Real > &v, const Vector< Real > &x, Real eps=0)
Set variables to zero if they correspond to the -active set.
ESecant StringToESecant(std::string s)
const ROL::Ptr< Vector< Real > > x_
ROL::Ptr< Secant< Real > > secant_
Secant object.
virtual void zero()
Set to zero vector.
Real scale_
Scale for dual variables in the active set, .
Defines the linear algebra or vector space interface.
ROL::Ptr< Vector< Real > > gtmp_
Container for temporary gradient storage.
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
const ROL::Ptr< Vector< Real > > xlam_
virtual void pruneLowerActive(Vector< Real > &v, const Vector< Real > &x, Real eps=0)
Set variables to zero if they correspond to the lower -active set.
ROL::Ptr< Vector< Real > > res_
Container for optimality system residual for quadratic model.
int flagCR_
CR termination flag.
int flag_
PDAS termination flag.
State for algorithm class. Will be used for restarts.
ROL::Ptr< Vector< Real > > x0_
Container for initial priaml variables.
const ROL::Ptr< Objective< Real > > obj_
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
PrimalDualActiveSetStep(ROL::ParameterList &parlist)
Constructor.
ESecant
Enumeration of secant update algorithms.
ROL::Ptr< Vector< Real > > rtmp_
Container for temporary right hand side storage.
ROL::Ptr< Vector< Real > > Ag_
Container for gradient projected onto active set.
ROL::Ptr< StepState< Real > > getState(void)
void initialize(Vector< Real > &x, const Vector< Real > &s, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Initialize step with bound constraint.
const ROL::Ptr< BoundConstraint< Real > > bnd_
bool feasible_
Flag whether the current iterate is feasible or not.
const ROL::Ptr< Objective< Real > > obj_
Provides interface for and implements limited-memory secant operators.
void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply linear operator.
Real computeCriticalityMeasure(Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con, Real tol)
Compute the gradient-based criticality measure.
virtual void pruneUpperActive(Vector< Real > &v, const Vector< Real > &x, Real eps=0)
Set variables to zero if they correspond to the upper -active set.
ROL::Ptr< Vector< Real > > iterateVec
ROL::Ptr< Vector< Real > > xbnd_
Container for primal variable bounds.
const ROL::Ptr< Secant< Real > > secant_
Provides the interface to apply a linear operator.
Provides the interface to apply upper and lower bound constraints.
void compute(Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Compute step.
Real gtol_
PDAS gradient stopping tolerance.
PrecondPD(const ROL::Ptr< Objective< Real > > &obj, const ROL::Ptr< BoundConstraint< Real > > &bnd, const ROL::Ptr< Vector< Real > > &x, const ROL::Ptr< Vector< Real > > &xlam, const Real eps=0, const ROL::Ptr< Secant< Real > > &secant=ROL::nullPtr, const bool useSecant=false)
Real stol_
PDAS minimum step size stopping tolerance.
ROL::Ptr< Vector< Real > > lambda_
Container for dual variables.
virtual void set(const Vector &x)
Set where .
virtual Real norm() const =0
Returns where .
ROL::Ptr< Krylov< Real > > krylov_
virtual const ROL::Ptr< const Vector< Real > > getLowerBound(void) const
Return the ref count pointer to the lower bound vector.
Real itol_
Inexact CR tolerance.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Real ROL_EPSILON(void)
Platform-dependent machine epsilon.
void applyInverse(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply inverse of linear operator.
ROL::Ptr< Vector< Real > > xlam_
Container for primal plus dual variables.
ESecant esec_
Enum for secant type.
const ROL::Ptr< Vector< Real > > x_
ROL::Ptr< Vector< Real > > v_
virtual bool isFeasible(const Vector< Real > &v)
Check if the vector, v, is feasible.
void update(Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Update step, if successful.
virtual std::string print(AlgorithmState< Real > &algo_state, bool print_header=false) const
Print iterate status.
ROL::Ptr< Vector< Real > > As_
Container for step projected onto active set.
ROL::Ptr< Vector< Real > > xtmp_
Container for temporary primal storage.
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
void apply(Vector< Real > &Hv, const Vector< Real > &v, Real &tol) const
Apply linear operator.
int maxit_
Maximum number of PDAS iterations.