ROL
ROL_AugmentedLagrangian_SimOpt.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_AUGMENTEDLAGRANGIAN_SIMOPT_H
45 #define ROL_AUGMENTEDLAGRANGIAN_SIMOPT_H
46 
47 #include "ROL_Objective_SimOpt.hpp"
50 #include "ROL_Vector.hpp"
51 #include "ROL_Types.hpp"
52 #include "ROL_Ptr.hpp"
53 #include <iostream>
54 
94 namespace ROL {
95 
96 template <class Real>
98 private:
99  // Required for Augmented Lagrangian definition
100  const ROL::Ptr<Objective_SimOpt<Real> > obj_;
101  ROL::Ptr<QuadraticPenalty_SimOpt<Real> > pen_;
103 
104  // Auxiliary storage
105  ROL::Ptr<Vector<Real> > dualSimVector_;
106  ROL::Ptr<Vector<Real> > dualOptVector_;
107 
108  // Objective and constraint evaluations
109  Real fval_;
110  ROL::Ptr<Vector<Real> > gradient1_;
111  ROL::Ptr<Vector<Real> > gradient2_;
112 
113  // Evaluation counters
114  int nfval_;
115  int ngval_;
116 
117  // User defined options
119 
120  // Flags to recompute quantities
124 
125 public:
127  const ROL::Ptr<Constraint_SimOpt<Real> > &con,
128  const Vector<Real> &multiplier,
129  const Real penaltyParameter,
130  const Vector<Real> &simVec,
131  const Vector<Real> &optVec,
132  const Vector<Real> &conVec,
133  ROL::ParameterList &parlist)
134  : obj_(obj), penaltyParameter_(penaltyParameter),
135  fval_(0), nfval_(0), ngval_(0), isValueComputed_(false),
137 
138  gradient1_ = simVec.dual().clone();
139  gradient2_ = optVec.dual().clone();
140  dualSimVector_ = simVec.dual().clone();
141  dualOptVector_ = optVec.dual().clone();
142 
143  ROL::ParameterList& sublist = parlist.sublist("Step").sublist("Augmented Lagrangian");
144  scaleLagrangian_ = sublist.get("Use Scaled Augmented Lagrangian", false);
145  int HessianApprox = sublist.get("Level of Hessian Approximation", 0);
146 
147  pen_ = ROL::makePtr<QuadraticPenalty_SimOpt<Real>>(con,multiplier,penaltyParameter,simVec,optVec,conVec,scaleLagrangian_,HessianApprox);
148  }
149 
150  virtual void update( const Vector<Real> &u, const Vector<Real> &z, bool flag = true, int iter = -1 ) {
151  obj_->update(u,z,flag,iter);
152  pen_->update(u,z,flag,iter);
153  isValueComputed_ = (flag ? false : isValueComputed_);
154  isGradient1Computed_ = (flag ? false : isGradient1Computed_);
155  isGradient2Computed_ = (flag ? false : isGradient2Computed_);
156  }
157 
158  virtual Real value( const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
159  // Compute objective function value
160  if ( !isValueComputed_ ) {
161  fval_ = obj_->value(u,z,tol); nfval_++;
162  isValueComputed_ = true;
163  }
164  // Compute penalty term
165  Real pval = pen_->value(u,z,tol);
166  // Compute augmented Lagrangian
167  Real val = fval_;
168  if (scaleLagrangian_) {
169  val /= penaltyParameter_;
170  }
171  return val + pval;
172  }
173 
174  virtual void gradient_1( Vector<Real> &g, const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
175  // Compute objective function gradient
176  if ( !isGradient1Computed_ ) {
177  obj_->gradient_1(*gradient1_,u,z,tol); ngval_++;
178  isGradient1Computed_ = true;
179  }
180  g.set(*gradient1_);
181  // Compute gradient of penalty
182  pen_->gradient_1(*dualSimVector_,u,z,tol);
183  // Compute gradient of Augmented Lagrangian
184  if ( scaleLagrangian_ ) {
185  g.scale(static_cast<Real>(1)/penaltyParameter_);
186  }
187  g.plus(*dualSimVector_);
188  }
189 
190  virtual void gradient_2( Vector<Real> &g, const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
191  // Compute objective function gradient
192  if ( !isGradient2Computed_ ) {
193  obj_->gradient_2(*gradient2_,u,z,tol); ngval_++;
194  isGradient2Computed_ = true;
195  }
196  g.set(*gradient2_);
197  // Compute gradient of penalty
198  pen_->gradient_2(*dualOptVector_,u,z,tol);
199  // Compute gradient of Augmented Lagrangian
200  if ( scaleLagrangian_ ) {
201  g.scale(static_cast<Real>(1)/penaltyParameter_);
202  }
203  g.plus(*dualOptVector_);
204  }
205 
206  virtual void hessVec_11( Vector<Real> &hv, const Vector<Real> &v,
207  const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
208  // Apply objective Hessian to a vector
209  obj_->hessVec_11(hv,v,u,z,tol);
210  // Apply penalty Hessian to a vector
211  pen_->hessVec_11(*dualSimVector_,v,u,z,tol);
212  // Build hessVec of Augmented Lagrangian
213  if ( scaleLagrangian_ ) {
214  hv.scale(static_cast<Real>(1)/penaltyParameter_);
215  }
216  hv.plus(*dualSimVector_);
217  }
218 
219  virtual void hessVec_12( Vector<Real> &hv, const Vector<Real> &v,
220  const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
221  // Apply objective Hessian to a vector
222  obj_->hessVec_12(hv,v,u,z,tol);
223  // Apply penalty Hessian to a vector
224  pen_->hessVec_12(*dualSimVector_,v,u,z,tol);
225  // Build hessVec of Augmented Lagrangian
226  if ( scaleLagrangian_ ) {
227  hv.scale(static_cast<Real>(1)/penaltyParameter_);
228  }
229  hv.plus(*dualSimVector_);
230  }
231 
232  virtual void hessVec_21( Vector<Real> &hv, const Vector<Real> &v,
233  const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
234  // Apply objective Hessian to a vector
235  obj_->hessVec_21(hv,v,u,z,tol);
236  // Apply penalty Hessian to a vector
237  pen_->hessVec_21(*dualOptVector_,v,u,z,tol);
238  // Build hessVec of Augmented Lagrangian
239  if ( scaleLagrangian_ ) {
240  hv.scale(static_cast<Real>(1)/penaltyParameter_);
241  }
242  hv.plus(*dualOptVector_);
243  }
244 
245  virtual void hessVec_22( Vector<Real> &hv, const Vector<Real> &v,
246  const Vector<Real> &u, const Vector<Real> &z, Real &tol ) {
247  // Apply objective Hessian to a vector
248  obj_->hessVec_22(hv,v,u,z,tol);
249  // Apply penalty Hessian to a vector
250  pen_->hessVec_22(*dualOptVector_,v,u,z,tol);
251  // Build hessVec of Augmented Lagrangian
252  if ( scaleLagrangian_ ) {
253  hv.scale(static_cast<Real>(1)/penaltyParameter_);
254  }
255  hv.plus(*dualOptVector_);
256  }
257 
258  // Return objective function value
259  virtual Real getObjectiveValue(const Vector<Real> &u, const Vector<Real> &z) {
260  Real tol = std::sqrt(ROL_EPSILON<Real>());
261  // Evaluate objective function value
262  if ( !isValueComputed_ ) {
263  fval_ = obj_->value(u,z,tol); nfval_++;
264  isValueComputed_ = true;
265  }
266  return fval_;
267  }
268 
269  // Return constraint value
270  virtual void getConstraintVec(Vector<Real> &c, const Vector<Real> &u, const Vector<Real> &z) {
271  pen_->getConstraintVec(c,u,z);
272  }
273 
274  // Return total number of constraint evaluations
275  virtual int getNumberConstraintEvaluations(void) const {
276  return pen_->getNumberConstraintEvaluations();
277  }
278 
279  // Return total number of objective evaluations
280  virtual int getNumberFunctionEvaluations(void) const {
281  return nfval_;
282  }
283 
284  // Return total number of gradient evaluations
285  virtual int getNumberGradientEvaluations(void) const {
286  return ngval_;
287  }
288 
289  // Reset with upated penalty parameter
290  virtual void reset(const Vector<Real> &multiplier, const Real penaltyParameter) {
291  nfval_ = 0; ngval_ = 0;
292  pen_->reset(multiplier,penaltyParameter);
293  }
294 }; // class AugmentedLagrangian
295 
296 } // namespace ROL
297 
298 #endif
Provides the interface to evaluate simulation-based objective functions.
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 hessVec_22(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
virtual void scale(const Real alpha)=0
Compute where .
virtual void plus(const Vector &x)=0
Compute , where .
AugmentedLagrangian_SimOpt(const ROL::Ptr< Objective_SimOpt< Real > > &obj, const ROL::Ptr< Constraint_SimOpt< Real > > &con, const Vector< Real > &multiplier, const Real penaltyParameter, const Vector< Real > &simVec, const Vector< Real > &optVec, const Vector< Real > &conVec, ROL::ParameterList &parlist)
virtual void hessVec_12(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Contains definitions of custom data types in ROL.
virtual void hessVec_21(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
virtual void update(const Vector< Real > &u, const Vector< Real > &z, bool flag=true, int iter=-1)
Update objective function. u is an iterate, z is an iterate, flag = true if the iterate has changed...
ROL::Ptr< QuadraticPenalty_SimOpt< Real > > pen_
virtual Real value(const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute value.
Provides the interface to evaluate the SimOpt augmented Lagrangian.
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
const ROL::Ptr< Objective_SimOpt< Real > > obj_
virtual void hessVec_11(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Apply Hessian approximation to vector.
virtual void gradient_2(Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute gradient with respect to second component.
virtual void gradient_1(Vector< Real > &g, const Vector< Real > &u, const Vector< Real > &z, Real &tol)
Compute gradient with respect to first component.
virtual void reset(const Vector< Real > &multiplier, const Real penaltyParameter)
virtual Real getObjectiveValue(const Vector< Real > &u, const Vector< Real > &z)
virtual void set(const Vector &x)
Set where .
Definition: ROL_Vector.hpp:209
Defines the constraint operator interface for simulation-based optimization.
virtual void getConstraintVec(Vector< Real > &c, const Vector< Real > &u, const Vector< Real > &z)