ROL
ROL_PH_RiskObjective.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 PH_RISKOBJECTIVE_H
45 #define PH_RISKOBJECTIVE_H
46 
47 #include "ROL_Objective.hpp"
49 
56 namespace ROL {
57 
58 template <class Real>
59 class PH_RiskObjective : public Objective<Real> {
60 private:
61  const Ptr<Objective<Real>> obj_;
62  Ptr<ExpectationQuad<Real>> quad_;
63 
65  Real val_;
66 
69  Ptr<Vector<Real>> g_;
70 
71  void getValue(const Vector<Real> &x, Real &tol) {
72  if (!isValueComputed_) {
73  val_ = obj_->value(x,tol);
74  isValueComputed_ = true;
75  }
76  }
77 
78  void getGradient(const Vector<Real> &x, Real &tol) {
80  g_ = x.dual().clone();
82  }
83  if (!isGradientComputed_) {
84  obj_->gradient(*g_,x,tol);
85  isGradientComputed_ = true;
86  }
87  }
88 
89  Ptr<const Vector<Real>> getConstVector(const Vector<Real> &x) const {
90  const RiskVector<Real> &xrv = dynamic_cast<const RiskVector<Real>&>(x);
91  return xrv.getVector();
92  }
93 
94  Ptr<Vector<Real>> getVector(Vector<Real> &x) const {
95  RiskVector<Real> &xrv = dynamic_cast<RiskVector<Real>&>(x);
96  return xrv.getVector();
97  }
98 
99  Ptr<const std::vector<Real>> getConstStat(const Vector<Real> &x) const {
100  const RiskVector<Real> &xrv = dynamic_cast<const RiskVector<Real>&>(x);
101  Ptr<const std::vector<Real>> xstat = xrv.getStatistic();
102  if (xstat == nullPtr) {
103  xstat = makePtr<const std::vector<Real>>(0);
104  }
105  return xstat;
106  }
107 
108  Ptr<std::vector<Real>> getStat(Vector<Real> &x) const {
109  RiskVector<Real> &xrv = dynamic_cast<RiskVector<Real>&>(x);
110  Ptr<std::vector<Real>> xstat = xrv.getStatistic();
111  if (xstat == nullPtr) {
112  xstat = makePtr<std::vector<Real>>(0);
113  }
114  return xstat;
115  }
116 
117 public:
118 
120  ParameterList &parlist)
121  : obj_(obj),
122  isValueComputed_(false),
123  isGradientInitialized_(false),
124  isGradientComputed_(false) {
125  std::string risk = parlist.sublist("SOL").sublist("Risk Measure").get("Name","CVaR");
127  switch(ed) {
128  case RISKMEASURE_CVAR:
129  quad_ = makePtr<QuantileQuadrangle<Real>>(parlist); break;
131  quad_ = makePtr<MoreauYosidaCVaR<Real>>(parlist); break;
133  quad_ = makePtr<GenMoreauYosidaCVaR<Real>>(parlist); break;
135  quad_ = makePtr<LogExponentialQuadrangle<Real>>(parlist); break;
137  quad_ = makePtr<MeanVarianceQuadrangle<Real>>(parlist); break;
139  quad_ = makePtr<TruncatedMeanQuadrangle<Real>>(parlist); break;
141  quad_ = makePtr<LogQuantileQuadrangle<Real>>(parlist); break;
143  quad_ = makePtr<SmoothedWorstCaseQuadrangle<Real>>(parlist); break;
144 // case RISKMEASURE_CHI2DIVERGENCE:
145 // return makePtr<Chi2Divergence<Real>>(parlist);
146 // case RISKMEASURE_KLDIVERGENCE:
147 // return makePtr<KLDivergence<Real>>(parlist);
148  default:
149  ROL_TEST_FOR_EXCEPTION(true,std::invalid_argument,
150  "Invalid risk measure type " << risk << "!");
151  }
152  }
153 
154  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
155  Ptr<const Vector<Real>> xvec = getConstVector(x);
156  obj_->update(*xvec,flag,iter);
157  isValueComputed_ = false;
158  isGradientComputed_ = false;
159  }
160 
161  Real value( const Vector<Real> &x, Real &tol ) {
162  Ptr<const Vector<Real>> xvec = getConstVector(x);
163  Ptr<const std::vector<Real>> xstat = getConstStat(x);
164  getValue(*xvec,tol);
165  Real reg = quad_->regret(val_-(*xstat)[0],0);
166  return (*xstat)[0] + reg;
167  }
168 
169  void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
170  Ptr<Vector<Real>> gvec = getVector(g);
171  Ptr<std::vector<Real>> gstat = getStat(g);
172  Ptr<const Vector<Real>> xvec = getConstVector(x);
173  Ptr<const std::vector<Real>> xstat = getConstStat(x);
174  getValue(*xvec,tol);
175  Real reg = quad_->regret(val_-(*xstat)[0],1);
176  getGradient(*xvec,tol);
177  gvec->set(*g_); gvec->scale(reg);
178  (*gstat)[0] = static_cast<Real>(1)-reg;
179  }
180 
181  void hessVec( Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x, Real &tol ) {
182  Ptr<Vector<Real>> hvec = getVector(hv);
183  Ptr<std::vector<Real>> hstat = getStat(hv);
184  Ptr<const Vector<Real>> vvec = getConstVector(v);
185  Ptr<const std::vector<Real>> vstat = getConstStat(v);
186  Ptr<const Vector<Real>> xvec = getConstVector(x);
187  Ptr<const std::vector<Real>> xstat = getConstStat(x);
188  getValue(*xvec,tol);
189  Real reg1 = quad_->regret(val_-(*xstat)[0],1);
190  Real reg2 = quad_->regret(val_-(*xstat)[0],2);
191  getGradient(*xvec,tol);
192  Real gv = vvec->dot(g_->dual());
193  obj_->hessVec(*hvec,*vvec,*xvec,tol);
194  hvec->scale(reg1); hvec->axpy(reg2*(gv-(*vstat)[0]),*g_);
195  (*hstat)[0] = reg2*((*vstat)[0]-gv);
196  }
197 
198  void setParameter(const std::vector<Real> &param) {
199  obj_->setParameter(param);
201  }
202 
203 };
204 
205 }
206 #endif
Provides the interface to evaluate objective functions.
ROL::Ptr< std::vector< Real > > getStatistic(const int comp=0, const int index=0)
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
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
Ptr< Vector< Real > > getVector(Vector< Real > &x) const
ROL::Ptr< const Vector< Real > > getVector(void) const
Ptr< const Vector< Real > > getConstVector(const Vector< Real > &x) const
ERiskMeasure StringToERiskMeasure(std::string s)
void getGradient(const Vector< Real > &x, Real &tol)
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
Ptr< Vector< Real > > g_
PH_RiskObjective(const Ptr< Objective< Real >> &obj, ParameterList &parlist)
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
Provides the interface for the progressive hedging risk objective.
Ptr< const std::vector< Real > > getConstStat(const Vector< Real > &x) const
Ptr< ExpectationQuad< Real > > quad_
virtual void setParameter(const std::vector< Real > &param)
void getValue(const Vector< Real > &x, Real &tol)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Ptr< std::vector< Real > > getStat(Vector< Real > &x) const
Real value(const Vector< Real > &x, Real &tol)
Compute value.
void setParameter(const std::vector< Real > &param)
const Ptr< Objective< Real > > obj_