ROL
ROL_AugmentedLagrangianStep.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_AUGMENTEDLAGRANGIANSTEP_H
45 #define ROL_AUGMENTEDLAGRANGIANSTEP_H
46 
48 #include "ROL_Types.hpp"
49 #include "ROL_Algorithm.hpp"
50 #include "ROL_LineSearchStep.hpp"
51 #include "ROL_TrustRegionStep.hpp"
52 #include "ROL_ParameterList.hpp"
53 
141 namespace ROL {
142 
143 template <class Real>
144 class AugmentedLagrangianStep : public Step<Real> {
145 private:
146  ROL::Ptr<Algorithm<Real> > algo_;
147  ROL::Ptr<Vector<Real> > x_;
148  ROL::Ptr<BoundConstraint<Real> > bnd_;
149 
150  ROL::ParameterList parlist_;
151  // Lagrange multiplier update
158  // Optimality tolerance update
163  // Feasibility tolerance update
168  // Subproblem information
169  bool print_;
170  int maxit_;
172  std::string subStep_;
176  // Scaling information
178  Real fscale_;
179  Real cscale_;
180  // Verbosity flag
182 
184  const Real mu, Objective<Real> &obj,
185  BoundConstraint<Real> &bnd) {
187  = dynamic_cast<AugmentedLagrangian<Real>&>(obj);
188  Real gnorm = 0., tol = std::sqrt(ROL_EPSILON<Real>());
189  augLag.gradient(g,x,tol);
190  if ( scaleLagrangian_ ) {
191  g.scale(mu);
192  }
193  // Compute norm of projected gradient
194  if (bnd.isActivated()) {
195  x_->set(x);
196  x_->axpy(static_cast<Real>(-1),g.dual());
197  bnd.project(*x_);
198  x_->axpy(static_cast<Real>(-1),x);
199  gnorm = x_->norm();
200  }
201  else {
202  gnorm = g.norm();
203  }
204  return gnorm;
205  }
206 
207 public:
208 
210  using Step<Real>::compute;
211  using Step<Real>::update;
212 
214 
215  AugmentedLagrangianStep(ROL::ParameterList &parlist)
216  : Step<Real>(), algo_(ROL::nullPtr),
217  x_(ROL::nullPtr), parlist_(parlist), subproblemIter_(0) {
218  Real one(1), p1(0.1), p9(0.9), ten(1.e1), oe8(1.e8), oem8(1.e-8);
219  ROL::ParameterList& sublist = parlist.sublist("Step").sublist("Augmented Lagrangian");
220  useDefaultInitPen_ = sublist.get("Use Default Initial Penalty Parameter",true);
221  Step<Real>::getState()->searchSize = sublist.get("Initial Penalty Parameter",ten);
222  // Multiplier update parameters
223  scaleLagrangian_ = sublist.get("Use Scaled Augmented Lagrangian", false);
224  minPenaltyLowerBound_ = sublist.get("Penalty Parameter Reciprocal Lower Bound", p1);
226  penaltyUpdate_ = sublist.get("Penalty Parameter Growth Factor", ten);
227  maxPenaltyParam_ = sublist.get("Maximum Penalty Parameter", oe8);
228  // Optimality tolerance update
229  optIncreaseExponent_ = sublist.get("Optimality Tolerance Update Exponent", one);
230  optDecreaseExponent_ = sublist.get("Optimality Tolerance Decrease Exponent", one);
231  optToleranceInitial_ = sublist.get("Initial Optimality Tolerance", one);
232  // Feasibility tolerance update
233  feasIncreaseExponent_ = sublist.get("Feasibility Tolerance Update Exponent", p1);
234  feasDecreaseExponent_ = sublist.get("Feasibility Tolerance Decrease Exponent", p9);
235  feasToleranceInitial_ = sublist.get("Initial Feasibility Tolerance", one);
236  // Subproblem information
237  print_ = sublist.get("Print Intermediate Optimization History", false);
238  maxit_ = sublist.get("Subproblem Iteration Limit", 1000);
239  subStep_ = sublist.get("Subproblem Step Type", "Trust Region");
240  parlist_.sublist("Step").set("Type",subStep_);
241  parlist_.sublist("Status Test").set("Iteration Limit",maxit_);
242  // Verbosity setting
243  verbosity_ = parlist.sublist("General").get("Print Verbosity", 0);
244  print_ = (verbosity_ > 0 ? true : print_);
245  // Outer iteration tolerances
246  outerFeasTolerance_ = parlist.sublist("Status Test").get("Constraint Tolerance", oem8);
247  outerOptTolerance_ = parlist.sublist("Status Test").get("Gradient Tolerance", oem8);
248  outerStepTolerance_ = parlist.sublist("Status Test").get("Step Tolerance", oem8);
249  // Scaling
250  useDefaultScaling_ = sublist.get("Use Default Problem Scaling", true);
251  fscale_ = sublist.get("Objective Scaling", 1.0);
252  cscale_ = sublist.get("Constraint Scaling", 1.0);
253  }
254 
259  AlgorithmState<Real> &algo_state ) {
260  bnd_ = ROL::makePtr<BoundConstraint<Real>>();
261  bnd_->deactivate();
262  initialize(x,g,l,c,obj,con,*bnd_,algo_state);
263  }
264 
269  AlgorithmState<Real> &algo_state ) {
270  Real one(1), TOL(1.e-2);
272  = dynamic_cast<AugmentedLagrangian<Real>&>(obj);
273  // Initialize step state
274  ROL::Ptr<StepState<Real> > state = Step<Real>::getState();
275  state->descentVec = x.clone();
276  state->gradientVec = g.clone();
277  state->constraintVec = c.clone();
278  // Initialize additional storage
279  x_ = x.clone();
280  // Initialize the algorithm state
281  algo_state.nfval = 0;
282  algo_state.ncval = 0;
283  algo_state.ngrad = 0;
284  // Project x onto the feasible set
285  if ( bnd.isActivated() ) {
286  bnd.project(x);
287  }
288  bnd.update(x,true,algo_state.iter);
289  // Update objective and constraint.
290  augLag.update(x,true,algo_state.iter);
291  if (useDefaultScaling_) {
292  fscale_ = one/std::max(one,augLag.getObjectiveGradient(x)->norm());
293  try {
294  Real tol = std::sqrt(ROL_EPSILON<Real>());
295  Ptr<Vector<Real>> ji = x.clone();
296  Real maxji(0), normji(0);
297  for (int i = 0; i < c.dimension(); ++i) {
298  con.applyAdjointJacobian(*ji,*c.basis(i),x,tol);
299  normji = ji->norm();
300  maxji = std::max(normji,maxji);
301  }
302  cscale_ = one/std::max(one,maxji);
303  }
304  catch (std::exception &e) {
305  cscale_ = one;
306  }
307  }
308  augLag.setScaling(fscale_,cscale_);
309  algo_state.value = augLag.getObjectiveValue(x);
310  algo_state.gnorm = computeGradient(*(state->gradientVec),x,state->searchSize,obj,bnd);
311  augLag.getConstraintVec(*(state->constraintVec),x);
312  algo_state.cnorm = (state->constraintVec)->norm();
313  if (useDefaultInitPen_) {
314  Step<Real>::getState()->searchSize
315  = std::max(static_cast<Real>(1e-8),std::min(static_cast<Real>(10)*
316  std::max(one,std::abs(fscale_*algo_state.value))
317  /std::max(one,std::pow(cscale_*algo_state.cnorm,2)),
318  static_cast<Real>(1e-2)*maxPenaltyParam_));
319  }
320  // Update evaluation counters
321  algo_state.ncval += augLag.getNumberConstraintEvaluations();
322  algo_state.nfval += augLag.getNumberFunctionEvaluations();
323  algo_state.ngrad += augLag.getNumberGradientEvaluations();
324  // Initialize intermediate stopping tolerances
325  minPenaltyReciprocal_ = std::min(one/state->searchSize,minPenaltyLowerBound_);
326  optTolerance_ = std::max(TOL*outerOptTolerance_,
328  optTolerance_ = std::min(optTolerance_,TOL*algo_state.gnorm);
329  feasTolerance_ = std::max(TOL*outerFeasTolerance_,
331  if (verbosity_ > 0) {
332  std::cout << std::endl;
333  std::cout << "Augmented Lagrangian Initialize" << std::endl;
334  std::cout << "Objective Scaling: " << fscale_ << std::endl;
335  std::cout << "Constraint Scaling: " << cscale_ << std::endl;
336  std::cout << std::endl;
337  }
338  }
339 
342  void compute( Vector<Real> &s, const Vector<Real> &x, const Vector<Real> &l,
343  Objective<Real> &obj, Constraint<Real> &con,
344  AlgorithmState<Real> &algo_state ) {
345  compute(s,x,l,obj,con,*bnd_,algo_state);
346  }
347 
350  void compute( Vector<Real> &s, const Vector<Real> &x, const Vector<Real> &l,
351  Objective<Real> &obj, Constraint<Real> &con,
352  BoundConstraint<Real> &bnd, AlgorithmState<Real> &algo_state ) {
353  Real one(1);
355  = dynamic_cast<AugmentedLagrangian<Real>&>(obj);
356  parlist_.sublist("Status Test").set("Gradient Tolerance",optTolerance_);
357  parlist_.sublist("Status Test").set("Step Tolerance",1.e-6*optTolerance_);
358  algo_ = ROL::makePtr<Algorithm<Real>>(subStep_,parlist_,false);
359  x_->set(x);
360  if ( bnd.isActivated() ) {
361  algo_->run(*x_,augLag,bnd,print_);
362  }
363  else {
364  algo_->run(*x_,augLag,print_);
365  }
366  s.set(*x_); s.axpy(-one,x);
367  subproblemIter_ = (algo_->getState())->iter;
368  }
369 
374  AlgorithmState<Real> &algo_state ) {
375  update(x,l,s,obj,con,*bnd_,algo_state);
376  }
377 
383  AlgorithmState<Real> &algo_state ) {
384  Real one(1), oem2(1.e-2);
386  = dynamic_cast<AugmentedLagrangian<Real>&>(obj);
387  ROL::Ptr<StepState<Real> > state = Step<Real>::getState();
388  state->SPiter = subproblemIter_;
389  // Update the step and store in state
390  x.plus(s);
391  algo_state.iterateVec->set(x);
392  state->descentVec->set(s);
393  algo_state.snorm = s.norm();
394  algo_state.iter++;
395  // Update objective function value
396  algo_state.value = augLag.getObjectiveValue(x);
397  // Update constraint value
398  augLag.getConstraintVec(*(state->constraintVec),x);
399  algo_state.cnorm = (state->constraintVec)->norm();
400  // Compute gradient of the augmented Lagrangian
401  algo_state.gnorm = computeGradient(*(state->gradientVec),x,state->searchSize,obj,bnd);
402  algo_state.gnorm /= std::min(fscale_,cscale_);
403  // Update evaluation counters
404  algo_state.nfval += augLag.getNumberFunctionEvaluations();
405  algo_state.ngrad += augLag.getNumberGradientEvaluations();
406  algo_state.ncval += augLag.getNumberConstraintEvaluations();
407  // Update objective function and constraints
408  augLag.update(x,true,algo_state.iter);
409  bnd.update(x,true,algo_state.iter);
410  // Update multipliers
411  minPenaltyReciprocal_ = std::min(one/state->searchSize,minPenaltyLowerBound_);
412  if ( cscale_*algo_state.cnorm < feasTolerance_ ) {
413  l.axpy(state->searchSize*cscale_,(state->constraintVec)->dual());
414  optTolerance_ = std::max(oem2*outerOptTolerance_,
416  feasTolerance_ = std::max(oem2*outerFeasTolerance_,
418  // Update Algorithm State
419  algo_state.snorm += state->searchSize*cscale_*algo_state.cnorm;
420  algo_state.lagmultVec->set(l);
421  }
422  else {
423  state->searchSize = std::min(penaltyUpdate_*state->searchSize,maxPenaltyParam_);
424  optTolerance_ = std::max(oem2*outerOptTolerance_,
426  feasTolerance_ = std::max(oem2*outerFeasTolerance_,
428  }
429  augLag.reset(l,state->searchSize);
430  }
431 
434  std::string printHeader( void ) const {
435  std::stringstream hist;
436 
437  if(verbosity_>0) {
438  hist << std::string(114,'-') << std::endl;
439  hist << "Augmented Lagrangian status output definitions" << std::endl << std::endl;
440  hist << " iter - Number of iterates (steps taken)" << std::endl;
441  hist << " fval - Objective function value" << std::endl;
442  hist << " cnorm - Norm of the constraint violation" << std::endl;
443  hist << " gLnorm - Norm of the gradient of the Lagrangian" << std::endl;
444  hist << " snorm - Norm of the step" << std::endl;
445  hist << " penalty - Penalty parameter" << std::endl;
446  hist << " feasTol - Feasibility tolerance" << std::endl;
447  hist << " optTol - Optimality tolerance" << std::endl;
448  hist << " #fval - Number of times the objective was computed" << std::endl;
449  hist << " #grad - Number of times the gradient was computed" << std::endl;
450  hist << " #cval - Number of times the constraint was computed" << std::endl;
451  hist << " subIter - Number of iterations to solve subproblem" << std::endl;
452  hist << std::string(114,'-') << std::endl;
453  }
454  hist << " ";
455  hist << std::setw(6) << std::left << "iter";
456  hist << std::setw(15) << std::left << "fval";
457  hist << std::setw(15) << std::left << "cnorm";
458  hist << std::setw(15) << std::left << "gLnorm";
459  hist << std::setw(15) << std::left << "snorm";
460  hist << std::setw(10) << std::left << "penalty";
461  hist << std::setw(10) << std::left << "feasTol";
462  hist << std::setw(10) << std::left << "optTol";
463  hist << std::setw(8) << std::left << "#fval";
464  hist << std::setw(8) << std::left << "#grad";
465  hist << std::setw(8) << std::left << "#cval";
466  hist << std::setw(8) << std::left << "subIter";
467  hist << std::endl;
468  return hist.str();
469  }
470 
473  std::string printName( void ) const {
474  std::stringstream hist;
475  hist << std::endl << " Augmented Lagrangian Solver";
476  hist << std::endl;
477  hist << "Subproblem Solver: " << subStep_ << std::endl;
478  return hist.str();
479  }
480 
483  std::string print( AlgorithmState<Real> &algo_state, bool pHeader = false ) const {
484  std::stringstream hist;
485  hist << std::scientific << std::setprecision(6);
486  if ( algo_state.iter == 0 ) {
487  hist << printName();
488  }
489  if ( pHeader ) {
490  hist << printHeader();
491  }
492  if ( algo_state.iter == 0 ) {
493  hist << " ";
494  hist << std::setw(6) << std::left << algo_state.iter;
495  hist << std::setw(15) << std::left << algo_state.value;
496  hist << std::setw(15) << std::left << algo_state.cnorm;
497  hist << std::setw(15) << std::left << algo_state.gnorm;
498  hist << std::setw(15) << std::left << " ";
499  hist << std::scientific << std::setprecision(2);
500  hist << std::setw(10) << std::left << Step<Real>::getStepState()->searchSize;
501  hist << std::setw(10) << std::left << std::max(feasTolerance_,outerFeasTolerance_);
502  hist << std::setw(10) << std::left << std::max(optTolerance_,outerOptTolerance_);
503  hist << std::endl;
504  }
505  else {
506  hist << " ";
507  hist << std::setw(6) << std::left << algo_state.iter;
508  hist << std::setw(15) << std::left << algo_state.value;
509  hist << std::setw(15) << std::left << algo_state.cnorm;
510  hist << std::setw(15) << std::left << algo_state.gnorm;
511  hist << std::setw(15) << std::left << algo_state.snorm;
512  hist << std::scientific << std::setprecision(2);
513  hist << std::setw(10) << std::left << Step<Real>::getStepState()->searchSize;
514  hist << std::setw(10) << std::left << feasTolerance_;
515  hist << std::setw(10) << std::left << optTolerance_;
516  hist << std::scientific << std::setprecision(6);
517  hist << std::setw(8) << std::left << algo_state.nfval;
518  hist << std::setw(8) << std::left << algo_state.ngrad;
519  hist << std::setw(8) << std::left << algo_state.ncval;
520  hist << std::setw(8) << std::left << subproblemIter_;
521  hist << std::endl;
522  }
523  return hist.str();
524  }
525 
531  AlgorithmState<Real> &algo_state ) {}
532 
538  AlgorithmState<Real> &algo_state ) {}
539 
540 }; // class AugmentedLagrangianStep
541 
542 } // namespace ROL
543 
544 #endif
Provides the interface to evaluate objective functions.
Provides the interface to evaluate the augmented Lagrangian.
void initialize(Vector< Real > &x, const Vector< Real > &g, Vector< Real > &l, const Vector< Real > &c, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state)
Initialize step with equality constraint.
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis...
Definition: ROL_Vector.hpp:226
virtual void scale(const Real alpha)=0
Compute where .
virtual ROL::Ptr< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual int dimension() const
Return dimension of the vector space.
Definition: ROL_Vector.hpp:196
virtual ROL::Ptr< Vector > basis(const int i) const
Return i-th basis vector.
Definition: ROL_Vector.hpp:182
std::string printHeader(void) const
Print iterate header.
virtual void plus(const Vector &x)=0
Compute , where .
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:153
bool isActivated(void) const
Check if bounds are on.
virtual int getNumberConstraintEvaluations(void) const
Provides the interface to compute optimization steps.
Definition: ROL_Step.hpp:69
AugmentedLagrangianStep(ROL::ParameterList &parlist)
Contains definitions of custom data types in ROL.
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update bounds.
virtual Real getObjectiveValue(const Vector< Real > &x)
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
Provides the interface to compute augmented Lagrangian steps.
std::string printName(void) const
Print step name.
Real computeGradient(Vector< Real > &g, const Vector< Real > &x, const Real mu, Objective< Real > &obj, BoundConstraint< Real > &bnd)
ROL::Ptr< BoundConstraint< Real > > bnd_
virtual void reset(const Vector< Real > &multiplier, const Real penaltyParameter)
State for algorithm class. Will be used for restarts.
Definition: ROL_Types.hpp:143
ROL::Ptr< Algorithm< Real > > algo_
ROL::Ptr< StepState< Real > > getState(void)
Definition: ROL_Step.hpp:74
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
ROL::Ptr< Vector< Real > > iterateVec
Definition: ROL_Types.hpp:157
virtual void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
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).
Provides the interface to apply upper and lower bound constraints.
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 and bound constraints.
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 .
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).
virtual int getNumberGradientEvaluations(void) const
std::string print(AlgorithmState< Real > &algo_state, bool pHeader=false) const
Print iterate status.
ROL::Ptr< Vector< Real > > lagmultVec
Definition: ROL_Types.hpp:158
virtual int getNumberFunctionEvaluations(void) const
void setScaling(const Real fscale, const Real cscale=1.0)
const Ptr< const Vector< Real > > getObjectiveGradient(const Vector< Real > &x)
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
virtual Real norm() const =0
Returns where .
virtual void getConstraintVec(Vector< Real > &c, const Vector< Real > &x)
void update(Vector< Real > &x, Vector< Real > &l, const Vector< Real > &s, Objective< Real > &obj, Constraint< Real > &con, AlgorithmState< Real > &algo_state)
Update step, if successful (equality constraint).
void update(Vector< Real > &x, const Vector< Real > &s, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Update step, for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring.
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, AlgorithmState< Real > &algo_state)
Compute step (equality constraint).
void compute(Vector< Real > &s, const Vector< Real > &x, Objective< Real > &obj, BoundConstraint< Real > &con, AlgorithmState< Real > &algo_state)
Compute step for bound constraints; here only to satisfy the interface requirements, does nothing, needs refactoring.