44 #ifndef ROL_MOREAUYOSIDAPENALTYSTEP_H
45 #define ROL_MOREAUYOSIDAPENALTYSTEP_H
58 #include "ROL_ParameterList.hpp"
125 class AugmentedLagrangianStep;
127 template <
class Real>
128 class MoreauYosidaPenaltyStep :
public Step<Real> {
133 ROL::Ptr<Vector<Real>>
x_;
134 ROL::Ptr<Vector<Real>>
g_;
135 ROL::Ptr<Vector<Real>>
l_;
136 ROL::Ptr<BoundConstraint<Real>>
bnd_;
157 Real zerotol = std::sqrt(ROL_EPSILON<Real>());
164 myPen.
gradient(*(state->gradientVec), x, zerotol);
166 state->gradientVec->plus(*
g_);
167 gLnorm_ = (state->gradientVec)->norm();
169 con.
value(*(state->constraintVec),x, zerotol);
170 algo_state.
cnorm = (state->constraintVec)->norm();
185 Real zerotol = std::sqrt(ROL_EPSILON<Real>());
191 myPen.
gradient(*(state->gradientVec), x, zerotol);
192 gLnorm_ = (state->gradientVec)->norm();
194 algo_state.
cnorm =
static_cast<Real
>(0);
212 x_(ROL::nullPtr),
g_(ROL::nullPtr),
l_(ROL::nullPtr),
216 Real ten(10), oem6(1.e-6), oem8(1.e-8);
217 ROL::ParameterList& steplist = parlist.sublist(
"Step").sublist(
"Moreau-Yosida Penalty");
219 tau_ = steplist.get(
"Penalty Parameter Growth Factor",ten);
221 print_ = steplist.sublist(
"Subproblem").get(
"Print History",
false);
223 Real gtol = steplist.sublist(
"Subproblem").get(
"Optimality Tolerance",oem8);
224 Real ctol = steplist.sublist(
"Subproblem").get(
"Feasibility Tolerance",oem8);
225 Real stol = oem6*std::min(gtol,ctol);
226 int maxit = steplist.sublist(
"Subproblem").get(
"Iteration Limit",1000);
227 parlist_.sublist(
"Status Test").set(
"Gradient Tolerance", gtol);
228 parlist_.sublist(
"Status Test").set(
"Constraint Tolerance", ctol);
229 parlist_.sublist(
"Status Test").set(
"Step Tolerance", stol);
230 parlist_.sublist(
"Status Test").set(
"Iteration Limit", maxit);
232 stepname_ = steplist.sublist(
"Subproblem").get(
"Step Type",
"Composite Step");
244 state->descentVec = x.
clone();
245 state->gradientVec = g.
clone();
246 state->constraintVec = c.
clone();
256 algo_state.
nfval = 0;
257 algo_state.
ncval = 0;
258 algo_state.
ngrad = 0;
269 state->descentVec = x.
clone();
270 state->gradientVec = g.
clone();
279 algo_state.
nfval = 0;
280 algo_state.
ncval = 0;
281 algo_state.
ngrad = 0;
284 bnd_ = ROL::makePtr<BoundConstraint<Real>>();
297 Ptr<Objective<Real>> penObj;
299 Ptr<Objective<Real>> raw_obj = makePtrFromRef(obj);
300 Ptr<Constraint<Real>> raw_con = makePtrFromRef(con);
302 penObj = makePtr<AugmentedLagrangian<Real>>(raw_obj,raw_con,l,one,x,*(state->constraintVec),
parlist_);
306 Ptr<Objective<Real>> raw_obj = makePtrFromRef(obj);
307 Ptr<Constraint<Real>> raw_con = makePtrFromRef(con);
309 penObj = makePtr<Fletcher<Real>>(raw_obj,raw_con,x,*(state->constraintVec),
parlist_);
313 penObj = makePtrFromRef(obj);
320 x_->set(x);
l_->set(l);
364 state->descentVec->set(s);
376 state->searchSize *=
tau_;
381 algo_state.
ncval += (
algo_->getState())->ncval;
395 state->descentVec->set(s);
405 state->searchSize *=
tau_;
417 std::stringstream hist;
419 hist << std::setw(6) << std::left <<
"iter";
420 hist << std::setw(15) << std::left <<
"fval";
422 hist << std::setw(15) << std::left <<
"cnorm";
424 hist << std::setw(15) << std::left <<
"gnorm";
425 hist << std::setw(15) << std::left <<
"ifeas";
426 hist << std::setw(15) << std::left <<
"snorm";
427 hist << std::setw(10) << std::left <<
"penalty";
428 hist << std::setw(8) << std::left <<
"#fval";
429 hist << std::setw(8) << std::left <<
"#grad";
431 hist << std::setw(8) << std::left <<
"#cval";
433 hist << std::setw(8) << std::left <<
"subIter";
441 std::stringstream hist;
442 hist <<
"\n" <<
" Moreau-Yosida Penalty solver";
450 std::stringstream hist;
451 hist << std::scientific << std::setprecision(6);
452 if ( algo_state.
iter == 0 ) {
458 if ( algo_state.
iter == 0 ) {
460 hist << std::setw(6) << std::left << algo_state.
iter;
461 hist << std::setw(15) << std::left << algo_state.
value;
463 hist << std::setw(15) << std::left << algo_state.
cnorm;
465 hist << std::setw(15) << std::left <<
gLnorm_;
467 hist << std::setw(15) << std::left <<
" ";
468 hist << std::scientific << std::setprecision(2);
469 hist << std::setw(10) << std::left << Step<Real>::getStepState()->searchSize;
474 hist << std::setw(6) << std::left << algo_state.
iter;
475 hist << std::setw(15) << std::left << algo_state.
value;
477 hist << std::setw(15) << std::left << algo_state.
cnorm;
479 hist << std::setw(15) << std::left <<
gLnorm_;
481 hist << std::setw(15) << std::left << algo_state.
snorm;
482 hist << std::scientific << std::setprecision(2);
483 hist << std::setw(10) << std::left << Step<Real>::getStepState()->searchSize;
484 hist << std::scientific << std::setprecision(6);
485 hist << std::setw(8) << std::left << algo_state.
nfval;
486 hist << std::setw(8) << std::left << algo_state.
ngrad;
488 hist << std::setw(8) << std::left << algo_state.
ncval;
Provides the interface to evaluate objective functions.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update constraint functions. x is the optimization variable, flag = true if optimization variable is ...
~MoreauYosidaPenaltyStep()
EStep StringToEStep(std::string s)
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
ROL::Ptr< Vector< Real > > g_
virtual void plus(const Vector &x)=0
Compute , where .
void updateState(const Vector< Real > &x, const Vector< Real > &l, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
bool isActivated(void) const
Check if bounds are on.
std::string printHeader(void) const
Print iterate header.
Provides the interface to compute optimization steps.
int getNumberGradientEvaluations(void)
MoreauYosidaPenaltyStep(ROL::ParameterList &parlist)
Contains definitions of custom data types in ROL.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
ROL::Ptr< Vector< Real > > x_
void initialize(Vector< Real > &x, const Vector< Real > &g, Vector< Real > &l, const Vector< Real > &c, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Initialize step with equality constraint.
void updateMultipliers(Real mu, const ROL::Vector< Real > &x)
Defines the linear algebra or vector space interface.
virtual void value(Vector< Real > &c, const Vector< Real > &x, Real &tol)=0
Evaluate the constraint operator at .
ROL::Ptr< Algorithm< Real > > algo_
ROL::Ptr< BoundConstraint< Real > > bnd_
State for algorithm class. Will be used for restarts.
ROL::Ptr< StatusTest< Real > > status_
void update(Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Update step, for bound constraints.
ROL::Ptr< StepState< Real > > getState(void)
ROL::ParameterList parlist_
Provides the interface to evaluate the Moreau-Yosida penalty function.
void update(Vector< Real > &x, Vector< Real > &l, const Vector< Real > &s, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Update step, if successful (equality and bound constraints).
ROL::Ptr< Vector< Real > > iterateVec
ROL::Ptr< Vector< Real > > l_
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
void compute(Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Compute step for bound constraints.
void initialize(Vector< Real > &x, const Vector< Real > &g, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Initialize step without equality constraint.
Provides the interface to apply upper and lower bound constraints.
int getNumberFunctionEvaluations(void)
virtual void applyAdjointJacobian(Vector< Real > &ajv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply the adjoint of the the constraint Jacobian at , , to vector .
std::string printName(void) const
Print step name.
ROL::Ptr< Vector< Real > > lagmultVec
void updateState(const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update Moreau-Yosida penalty function.
virtual void set(const Vector &x)
Set where .
virtual Real norm() const =0
Returns where .
ROL::Ptr< Step< Real > > step_
std::string print(AlgorithmState< Real > &algo_state, bool pHeader=false) const
Print iterate status.
EStep
Enumeration of step types.
Real testComplementarity(const ROL::Vector< Real > &x)
Defines the general constraint operator interface.
virtual void project(Vector< Real > &x)
Project optimization variables onto the bounds.
void compute(Vector< Real > &s, const Vector< Real > &x, const Vector< Real > &l, Objective< Real > &obj, Constraint< Real > &con, BoundConstraint< Real > &bnd, AlgorithmState< Real > &algo_state)
Compute step (equality and bound constraints).