ROL
ROL_Reduced_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_REDUCED_AUGMENTEDLAGRANGIAN_SIMOPT_H
45 #define ROL_REDUCED_AUGMENTEDLAGRANGIAN_SIMOPT_H
46 
47 #include "ROL_Objective_SimOpt.hpp"
52 #include "ROL_Vector.hpp"
53 #include "ROL_Types.hpp"
54 #include "ROL_Ptr.hpp"
55 #include <iostream>
56 
95 namespace ROL {
96 
97 template <class Real>
99 private:
100  ROL::Ptr<AugmentedLagrangian_SimOpt<Real> > augLagSimOpt_;
101  ROL::Ptr<Reduced_Objective_SimOpt<Real> > rAugLagSimOpt_;
102  ROL::Ptr<Vector<Real> > state_;
103 
104  // Evaluation counters
105  int ngval_;
106 
107 public:
109  const ROL::Ptr<Constraint_SimOpt<Real> > &redCon,
110  const ROL::Ptr<Constraint_SimOpt<Real> > &augCon,
111  const ROL::Ptr<Vector<Real> > &state,
112  const ROL::Ptr<Vector<Real> > &control,
113  const ROL::Ptr<Vector<Real> > &adjoint,
114  const ROL::Ptr<Vector<Real> > &augConVec,
115  const ROL::Ptr<Vector<Real> > &multiplier,
116  const Real penaltyParameter,
117  ROL::ParameterList &parlist) : state_(state),
118  ngval_(0) {
119 
120  augLagSimOpt_ = ROL::makePtr<AugmentedLagrangian_SimOpt<Real>>(obj,
121  augCon,
122  *multiplier,
123  penaltyParameter,
124  *state,
125  *control,
126  *augConVec,
127  parlist);
128  rAugLagSimOpt_ = ROL::makePtr<Reduced_Objective_SimOpt<Real>>(augLagSimOpt_,redCon,state,control,adjoint);
129  rAugLagSimOpt_->update(*control);
130  Real tol = 1e-8;
131  rAugLagSimOpt_->value(*control,tol);
132  }
133 
134  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
135  rAugLagSimOpt_->update(x,flag,iter);
136  }
137 
138  Real value( const Vector<Real> &x, Real &tol ) {
139  return rAugLagSimOpt_->value(x,tol);
140  }
141 
142  void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
143  ngval_++;
144  rAugLagSimOpt_->gradient(g,x,tol);
145  }
146 
147  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
148  rAugLagSimOpt_->hessVec(hv,v,x,tol);
149  }
150 
151  // Return objective function value
153  return augLagSimOpt_->getObjectiveValue(*state_,x);
154  }
155 
156  // Return constraint value
158  augLagSimOpt_->getConstraintVec(c,*state_,x);
159  }
160 
161  // Return total number of constraint evaluations
163  return augLagSimOpt_->getNumberConstraintEvaluations();
164  }
165 
166  // Return total number of objective evaluations
168  return augLagSimOpt_->getNumberFunctionEvaluations();
169  }
170 
171  // Return total number of gradient evaluations
173  return ngval_;
174  //return augLagSimOpt_->getNumberGradientEvaluations();
175  }
176 
177  // Reset with upated penalty parameter
178  void reset(const Vector<Real> &multiplier, const Real penaltyParameter) {
179  ngval_ = 0;
180  augLagSimOpt_->reset(multiplier,penaltyParameter);
181  }
182 
183 // For parametrized (stochastic) objective functions and constraints
184 public:
185  void setParameter(const std::vector<Real> &param) {
187  rAugLagSimOpt_->setParameter(param);
188  }
189 }; // class AugmentedLagrangian
190 
191 } // namespace ROL
192 
193 #endif
Provides the interface to evaluate simulation-based objective functions.
Provides the interface to evaluate the augmented Lagrangian.
void getConstraintVec(Vector< Real > &c, const Vector< Real > &x)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Contains definitions of custom data types in ROL.
ROL::Ptr< Reduced_Objective_SimOpt< Real > > rAugLagSimOpt_
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
virtual void setParameter(const std::vector< Real > &param)
ROL::Ptr< AugmentedLagrangian_SimOpt< Real > > augLagSimOpt_
Defines the constraint operator interface for simulation-based optimization.
Reduced_AugmentedLagrangian_SimOpt(const ROL::Ptr< Objective_SimOpt< Real > > &obj, const ROL::Ptr< Constraint_SimOpt< Real > > &redCon, const ROL::Ptr< Constraint_SimOpt< Real > > &augCon, const ROL::Ptr< Vector< Real > > &state, const ROL::Ptr< Vector< Real > > &control, const ROL::Ptr< Vector< Real > > &adjoint, const ROL::Ptr< Vector< Real > > &augConVec, const ROL::Ptr< Vector< Real > > &multiplier, const Real penaltyParameter, ROL::ParameterList &parlist)
Provides the interface to evaluate the reduced SimOpt augmented Lagrangian.
Real value(const Vector< Real > &x, Real &tol)
Compute value.
void reset(const Vector< Real > &multiplier, const Real penaltyParameter)