ROL
ROL_SpectralRisk.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_SPECTRALRISK_HPP
45 #define ROL_SPECTRALRISK_HPP
46 
47 #include "ROL_MixedCVaR.hpp"
49 
86 namespace ROL {
87 
88 template<class Real>
89 class SpectralRisk : public RandVarFunctional<Real> {
90 private:
91  ROL::Ptr<MixedCVaR<Real> > mqq_;
92  ROL::Ptr<PlusFunction<Real> > plusFunction_;
93 
94  std::vector<Real> wts_;
95  std::vector<Real> pts_;
96 
97  void checkInputs(ROL::Ptr<Distribution<Real> > &dist = ROL::nullPtr) const {
98  ROL_TEST_FOR_EXCEPTION(plusFunction_ == ROL::nullPtr, std::invalid_argument,
99  ">>> ERROR (ROL::SpectralRisk): PlusFunction pointer is null!");
100  if ( dist != ROL::nullPtr) {
101  Real lb = dist->lowerBound();
102  Real ub = dist->upperBound();
103  ROL_TEST_FOR_EXCEPTION(lb < static_cast<Real>(0), std::invalid_argument,
104  ">>> ERROR (ROL::SpectralRisk): Distribution lower bound less than zero!");
105  ROL_TEST_FOR_EXCEPTION(ub > static_cast<Real>(1), std::invalid_argument,
106  ">>> ERROR (ROL::SpectralRisk): Distribution upper bound greater than one!");
107  }
108  }
109 
110 protected:
111  void buildMixedQuantile(const std::vector<Real> &pts, const std::vector<Real> &wts,
112  const ROL::Ptr<PlusFunction<Real> > &pf) {
113  pts_.clear(); pts_.assign(pts.begin(),pts.end());
114  wts_.clear(); wts_.assign(wts.begin(),wts.end());
115  plusFunction_ = pf;
116  mqq_ = ROL::makePtr<MixedCVaR<Real>>(pts,wts,pf);
117  }
118 
119  void buildQuadFromDist(std::vector<Real> &pts, std::vector<Real> &wts,
120  const int nQuad, const ROL::Ptr<Distribution<Real> > &dist) const {
121  const Real lo = dist->lowerBound(), hi = dist->upperBound();
122  const Real half(0.5), one(1), N(nQuad);
123  wts.clear(); wts.resize(nQuad);
124  pts.clear(); pts.resize(nQuad);
125  if ( hi >= one ) {
126  wts[0] = half/(N-half);
127  pts[0] = lo;
128  for (int i = 1; i < nQuad; ++i) {
129  wts[i] = one/(N-half);
130  pts[i] = dist->invertCDF(static_cast<Real>(i)/N);
131  }
132  }
133  else {
134  wts[0] = half/(N-one);
135  pts[0] = lo;
136  for (int i = 1; i < nQuad-1; ++i) {
137  wts[i] = one/(N-one);
138  pts[i] = dist->invertCDF(static_cast<Real>(i)/N);
139  }
140  wts[nQuad-1] = half/(N-one);
141  pts[nQuad-1] = hi;
142  }
143  }
144 
145  void printQuad(const std::vector<Real> &pts,
146  const std::vector<Real> &wts,
147  const bool print = false) const {
148  if ( print ) {
149  const int nQuad = wts.size();
150  std::cout << std::endl;
151  std::cout << std::scientific << std::setprecision(15);
152  std::cout << std::setw(25) << std::left << "Points"
153  << std::setw(25) << std::left << "Weights"
154  << std::endl;
155  for (int i = 0; i < nQuad; ++i) {
156  std::cout << std::setw(25) << std::left << pts[i]
157  << std::setw(25) << std::left << wts[i]
158  << std::endl;
159  }
160  std::cout << std::endl;
161  }
162  }
163 
164 
165 public:
166  SpectralRisk(void) : RandVarFunctional<Real>() {}
167 
168  SpectralRisk( const ROL::Ptr<Distribution<Real> > &dist,
169  const int nQuad,
170  const ROL::Ptr<PlusFunction<Real> > &pf)
171  : RandVarFunctional<Real>() {
172  // Build generalized trapezoidal rule
173  std::vector<Real> wts(nQuad), pts(nQuad);
174  buildQuadFromDist(pts,wts,nQuad,dist);
175  // Build mixed quantile quadrangle risk measure
176  buildMixedQuantile(pts,wts,pf);
177  // Check inputs
178  checkInputs(dist);
179  }
180 
181  SpectralRisk(ROL::ParameterList &parlist)
182  : RandVarFunctional<Real>() {
183  // Parse parameter list
184  ROL::ParameterList &list
185  = parlist.sublist("SOL").sublist("Risk Measure").sublist("Spectral Risk");
186  int nQuad = list.get("Number of Quadrature Points",5);
187  bool print = list.get("Print Quadrature to Screen",false);
188  // Build distribution
189  ROL::Ptr<Distribution<Real> > dist = DistributionFactory<Real>(list);
190  // Build plus function approximation
191  ROL::Ptr<PlusFunction<Real> > pf = ROL::makePtr<PlusFunction<Real>>(list);
192  // Build generalized trapezoidal rule
193  std::vector<Real> wts(nQuad), pts(nQuad);
194  buildQuadFromDist(pts,wts,nQuad,dist);
195  printQuad(pts,wts,print);
196  // Build mixed quantile quadrangle risk measure
197  buildMixedQuantile(pts,wts,pf);
198  // Check inputs
199  checkInputs(dist);
200  }
201 
202  SpectralRisk( const std::vector<Real> &pts, const std::vector<Real> &wts,
203  const ROL::Ptr<PlusFunction<Real> > &pf)
204  : RandVarFunctional<Real>() {
205  buildMixedQuantile(pts,wts,pf);
206  // Check inputs
207  checkInputs();
208  }
209 
210  Real computeStatistic(const std::vector<Real> &xstat) {
211  Real stat(0);
212  int nQuad = static_cast<int>(wts_.size());
213  for (int i = 0; i < nQuad; ++i) {
214  stat += wts_[i] * xstat[i];
215  }
216  return stat;
217  }
218 
219  void setStorage(const Ptr<SampledScalar<Real>> &value_storage,
220  const Ptr<SampledVector<Real>> &gradient_storage) {
221  RandVarFunctional<Real>::setStorage(value_storage,gradient_storage);
222  mqq_->setStorage(value_storage,gradient_storage);
223  }
224 
225  void setHessVecStorage(const Ptr<SampledScalar<Real>> &gradvec_storage,
226  const Ptr<SampledVector<Real>> &hessvec_storage) {
227  RandVarFunctional<Real>::setHessVecStorage(gradvec_storage,hessvec_storage);
228  mqq_->setHessVecStorage(gradvec_storage,hessvec_storage);
229  }
230 
231  void setSample(const std::vector<Real> &point, const Real weight) {
233  mqq_->setSample(point,weight);
234  }
235 
236  void resetStorage(bool flag = true) {
238  mqq_->resetStorage(flag);
239  }
240 
241  void initialize(const Vector<Real> &x) {
243  mqq_->initialize(x);
244  }
245 
246  Real computeStatistic(const Ptr<std::vector<Real>> &xstat) const {
247  return mqq_->computeStatistic(xstat);
248  }
249 
251  const Vector<Real> &x,
252  const std::vector<Real> &xstat,
253  Real &tol) {
254  mqq_->updateValue(obj,x,xstat,tol);
255  }
256 
258  const Vector<Real> &x,
259  const std::vector<Real> &xstat,
260  Real &tol) {
261  mqq_->updateGradient(obj,x,xstat,tol);
262  }
263 
265  const Vector<Real> &v,
266  const std::vector<Real> &vstat,
267  const Vector<Real> &x,
268  const std::vector<Real> &xstat,
269  Real &tol) {
270  mqq_->updateHessVec(obj,v,vstat,x,xstat,tol);
271  }
272 
273  Real getValue(const Vector<Real> &x,
274  const std::vector<Real> &xstat,
275  SampleGenerator<Real> &sampler) {
276  return mqq_->getValue(x,xstat,sampler);
277  }
278 
280  std::vector<Real> &gstat,
281  const Vector<Real> &x,
282  const std::vector<Real> &xstat,
283  SampleGenerator<Real> &sampler) {
284  mqq_->getGradient(g,gstat,x,xstat,sampler);
285  }
286 
288  std::vector<Real> &hvstat,
289  const Vector<Real> &v,
290  const std::vector<Real> &vstat,
291  const Vector<Real> &x,
292  const std::vector<Real> &xstat,
293  SampleGenerator<Real> &sampler) {
294  mqq_->getHessVec(hv,hvstat,v,vstat,x,xstat,sampler);
295  }
296 };
297 
298 }
299 
300 #endif
virtual void setHessVecStorage(const Ptr< SampledScalar< Real >> &gradvec_storage, const Ptr< SampledVector< Real >> &hessvec_storage)
Provides the interface to evaluate objective functions.
void setStorage(const Ptr< SampledScalar< Real >> &value_storage, const Ptr< SampledVector< Real >> &gradient_storage)
virtual void setSample(const std::vector< Real > &point, const Real weight)
void updateHessVec(Objective< Real > &obj, const Vector< Real > &v, const std::vector< Real > &vstat, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal risk measure storage for Hessian-time-a-vector computation.
void buildMixedQuantile(const std::vector< Real > &pts, const std::vector< Real > &wts, const ROL::Ptr< PlusFunction< Real > > &pf)
std::vector< Real > pts_
void checkInputs(ROL::Ptr< Distribution< Real > > &dist=ROL::nullPtr) const
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
void getHessVec(Vector< Real > &hv, std::vector< Real > &hvstat, const Vector< Real > &v, const std::vector< Real > &vstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure Hessian-times-a-vector.
virtual void resetStorage(bool flag=true)
Reset internal storage.
virtual void setStorage(const Ptr< SampledScalar< Real >> &value_storage, const Ptr< SampledVector< Real >> &gradient_storage)
void resetStorage(bool flag=true)
Reset internal storage.
void updateGradient(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal risk measure storage for gradient computation.
void setSample(const std::vector< Real > &point, const Real weight)
ROL::Ptr< MixedCVaR< Real > > mqq_
Provides an interface for spectral risk measures.
Real computeStatistic(const std::vector< Real > &xstat)
void buildQuadFromDist(std::vector< Real > &pts, std::vector< Real > &wts, const int nQuad, const ROL::Ptr< Distribution< Real > > &dist) const
ROL::Ptr< PlusFunction< Real > > plusFunction_
SpectralRisk(ROL::ParameterList &parlist)
void getGradient(Vector< Real > &g, std::vector< Real > &gstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure (sub)gradient.
SpectralRisk(const ROL::Ptr< Distribution< Real > > &dist, const int nQuad, const ROL::Ptr< PlusFunction< Real > > &pf)
void updateValue(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal storage for value computation.
Real getValue(const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure value.
void printQuad(const std::vector< Real > &pts, const std::vector< Real > &wts, const bool print=false) const
std::vector< Real > wts_
void setHessVecStorage(const Ptr< SampledScalar< Real >> &gradvec_storage, const Ptr< SampledVector< Real >> &hessvec_storage)
Real computeStatistic(const Ptr< std::vector< Real >> &xstat) const
Provides the interface to implement any functional that maps a random variable to a (extended) real n...
void initialize(const Vector< Real > &x)
Initialize temporary variables.
SpectralRisk(const std::vector< Real > &pts, const std::vector< Real > &wts, const ROL::Ptr< PlusFunction< Real > > &pf)
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.