ROL
ROL_SpectralRisk.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Rapid Optimization Library (ROL) Package
4//
5// Copyright 2014 NTESS and the ROL contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef ROL_SPECTRALRISK_HPP
11#define ROL_SPECTRALRISK_HPP
12
13#include "ROL_MixedCVaR.hpp"
15
52namespace ROL {
53
54template<class Real>
55class SpectralRisk : public RandVarFunctional<Real> {
56private:
57 ROL::Ptr<MixedCVaR<Real> > mqq_;
58 ROL::Ptr<PlusFunction<Real> > plusFunction_;
59
60 std::vector<Real> wts_;
61 std::vector<Real> pts_;
62
63 void checkInputs(const ROL::Ptr<Distribution<Real> > &dist = ROL::nullPtr) const {
64 ROL_TEST_FOR_EXCEPTION(plusFunction_ == ROL::nullPtr, std::invalid_argument,
65 ">>> ERROR (ROL::SpectralRisk): PlusFunction pointer is null!");
66 if ( dist != ROL::nullPtr) {
67 Real lb = dist->lowerBound();
68 Real ub = dist->upperBound();
69 ROL_TEST_FOR_EXCEPTION(lb < static_cast<Real>(0), std::invalid_argument,
70 ">>> ERROR (ROL::SpectralRisk): Distribution lower bound less than zero!");
71 ROL_TEST_FOR_EXCEPTION(ub > static_cast<Real>(1), std::invalid_argument,
72 ">>> ERROR (ROL::SpectralRisk): Distribution upper bound greater than one!");
73 }
74 }
75
76protected:
77 void buildMixedQuantile(const std::vector<Real> &pts, const std::vector<Real> &wts,
78 const ROL::Ptr<PlusFunction<Real> > &pf) {
79 pts_.clear(); pts_.assign(pts.begin(),pts.end());
80 wts_.clear(); wts_.assign(wts.begin(),wts.end());
81 plusFunction_ = pf;
82 mqq_ = ROL::makePtr<MixedCVaR<Real>>(pts,wts,pf);
83 }
84
85 void buildQuadFromDist(std::vector<Real> &pts, std::vector<Real> &wts,
86 const int nQuad, const ROL::Ptr<Distribution<Real> > &dist) const {
87 const Real lo = dist->lowerBound(), hi = dist->upperBound();
88 const Real half(0.5), one(1), N(nQuad);
89 wts.clear(); wts.resize(nQuad);
90 pts.clear(); pts.resize(nQuad);
91 if ( hi >= one ) {
92 wts[0] = half/(N-half);
93 pts[0] = lo;
94 for (int i = 1; i < nQuad; ++i) {
95 wts[i] = one/(N-half);
96 pts[i] = dist->invertCDF(static_cast<Real>(i)/N);
97 }
98 }
99 else {
100 wts[0] = half/(N-one);
101 pts[0] = lo;
102 for (int i = 1; i < nQuad-1; ++i) {
103 wts[i] = one/(N-one);
104 pts[i] = dist->invertCDF(static_cast<Real>(i)/N);
105 }
106 wts[nQuad-1] = half/(N-one);
107 pts[nQuad-1] = hi;
108 }
109 }
110
111 void printQuad(const std::vector<Real> &pts,
112 const std::vector<Real> &wts,
113 const bool print = false) const {
114 if ( print ) {
115 const int nQuad = wts.size();
116 std::cout << std::endl;
117 std::cout << std::scientific << std::setprecision(15);
118 std::cout << std::setw(25) << std::left << "Points"
119 << std::setw(25) << std::left << "Weights"
120 << std::endl;
121 for (int i = 0; i < nQuad; ++i) {
122 std::cout << std::setw(25) << std::left << pts[i]
123 << std::setw(25) << std::left << wts[i]
124 << std::endl;
125 }
126 std::cout << std::endl;
127 }
128 }
129
130
131public:
133
134 SpectralRisk( const ROL::Ptr<Distribution<Real> > &dist,
135 const int nQuad,
136 const ROL::Ptr<PlusFunction<Real> > &pf)
137 : RandVarFunctional<Real>() {
138 // Build generalized trapezoidal rule
139 std::vector<Real> wts(nQuad), pts(nQuad);
140 buildQuadFromDist(pts,wts,nQuad,dist);
141 // Build mixed quantile quadrangle risk measure
142 buildMixedQuantile(pts,wts,pf);
143 // Check inputs
144 checkInputs(dist);
145 }
146
147 SpectralRisk(ROL::ParameterList &parlist)
148 : RandVarFunctional<Real>() {
149 // Parse parameter list
150 ROL::ParameterList &list
151 = parlist.sublist("SOL").sublist("Risk Measure").sublist("Spectral Risk");
152 int nQuad = list.get("Number of Quadrature Points",5);
153 bool print = list.get("Print Quadrature to Screen",false);
154 // Build distribution
155 ROL::Ptr<Distribution<Real> > dist = DistributionFactory<Real>(list);
156 // Build plus function approximation
157 ROL::Ptr<PlusFunction<Real> > pf = ROL::makePtr<PlusFunction<Real>>(list);
158 // Build generalized trapezoidal rule
159 std::vector<Real> wts(nQuad), pts(nQuad);
160 buildQuadFromDist(pts,wts,nQuad,dist);
161 printQuad(pts,wts,print);
162 // Build mixed quantile quadrangle risk measure
163 buildMixedQuantile(pts,wts,pf);
164 // Check inputs
165 checkInputs(dist);
166 }
167
168 SpectralRisk( const std::vector<Real> &pts, const std::vector<Real> &wts,
169 const ROL::Ptr<PlusFunction<Real> > &pf)
170 : RandVarFunctional<Real>() {
171 buildMixedQuantile(pts,wts,pf);
172 // Check inputs
173 checkInputs();
174 }
175
176 void setStorage(const Ptr<ScalarController<Real>> &value_storage,
177 const Ptr<VectorController<Real>> &gradient_storage) override {
178 RandVarFunctional<Real>::setStorage(value_storage,gradient_storage);
179 mqq_->setStorage(value_storage,gradient_storage);
180 }
181
182 void setHessVecStorage(const Ptr<ScalarController<Real>> &gradvec_storage,
183 const Ptr<VectorController<Real>> &hessvec_storage) override {
184 RandVarFunctional<Real>::setHessVecStorage(gradvec_storage,hessvec_storage);
185 mqq_->setHessVecStorage(gradvec_storage,hessvec_storage);
186 }
187
188 void setSample(const std::vector<Real> &point, const Real weight) override {
190 mqq_->setSample(point,weight);
191 }
192
193 void resetStorage(bool flag = true) override {
195 mqq_->resetStorage(flag);
196 }
197
198 void resetStorage(UpdateType type) override {
200 mqq_->resetStorage(type);
201 }
202
203 void initialize(const Vector<Real> &x) override {
205 mqq_->initialize(x);
206 }
207
208 Real computeStatistic(const Ptr<const std::vector<Real>> &xstat) const override {
209 return mqq_->computeStatistic(xstat);
210 }
211
213 const Vector<Real> &x,
214 const std::vector<Real> &xstat,
215 Real &tol) override {
216 mqq_->updateValue(obj,x,xstat,tol);
217 }
218
220 const Vector<Real> &x,
221 const std::vector<Real> &xstat,
222 Real &tol) override {
223 mqq_->updateGradient(obj,x,xstat,tol);
224 }
225
227 const Vector<Real> &v,
228 const std::vector<Real> &vstat,
229 const Vector<Real> &x,
230 const std::vector<Real> &xstat,
231 Real &tol) override {
232 mqq_->updateHessVec(obj,v,vstat,x,xstat,tol);
233 }
234
235 Real getValue(const Vector<Real> &x,
236 const std::vector<Real> &xstat,
237 SampleGenerator<Real> &sampler) override {
238 return mqq_->getValue(x,xstat,sampler);
239 }
240
242 std::vector<Real> &gstat,
243 const Vector<Real> &x,
244 const std::vector<Real> &xstat,
245 SampleGenerator<Real> &sampler) override {
246 mqq_->getGradient(g,gstat,x,xstat,sampler);
247 }
248
250 std::vector<Real> &hvstat,
251 const Vector<Real> &v,
252 const std::vector<Real> &vstat,
253 const Vector<Real> &x,
254 const std::vector<Real> &xstat,
255 SampleGenerator<Real> &sampler) override {
256 mqq_->getHessVec(hv,hvstat,v,vstat,x,xstat,sampler);
257 }
258};
259
260}
261
262#endif
Provides the interface to evaluate objective functions.
Provides the interface to implement any functional that maps a random variable to a (extended) real n...
virtual void setSample(const std::vector< Real > &point, const Real weight)
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.
virtual void setStorage(const Ptr< ScalarController< Real > > &value_storage, const Ptr< VectorController< Real > > &gradient_storage)
virtual void setHessVecStorage(const Ptr< ScalarController< Real > > &gradvec_storage, const Ptr< VectorController< Real > > &hessvec_storage)
virtual void resetStorage(bool flag=true)
Reset internal storage.
Provides an interface for spectral risk measures.
SpectralRisk(const ROL::Ptr< Distribution< Real > > &dist, const int nQuad, const ROL::Ptr< PlusFunction< Real > > &pf)
void printQuad(const std::vector< Real > &pts, const std::vector< Real > &wts, const bool print=false) const
std::vector< Real > wts_
void checkInputs(const ROL::Ptr< Distribution< Real > > &dist=ROL::nullPtr) const
void resetStorage(bool flag=true) override
Reset internal storage.
void buildQuadFromDist(std::vector< Real > &pts, std::vector< Real > &wts, const int nQuad, const ROL::Ptr< Distribution< Real > > &dist) const
Real computeStatistic(const Ptr< const std::vector< Real > > &xstat) const override
Compute statistic.
ROL::Ptr< PlusFunction< Real > > plusFunction_
void buildMixedQuantile(const std::vector< Real > &pts, const std::vector< Real > &wts, const ROL::Ptr< PlusFunction< Real > > &pf)
void initialize(const Vector< Real > &x) override
Initialize temporary variables.
void setStorage(const Ptr< ScalarController< Real > > &value_storage, const Ptr< VectorController< Real > > &gradient_storage) override
void updateValue(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol) override
Update internal storage for value computation.
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) override
Update internal risk measure storage for Hessian-time-a-vector computation.
void setSample(const std::vector< Real > &point, const Real weight) override
void getGradient(Vector< Real > &g, std::vector< Real > &gstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler) override
Return risk measure (sub)gradient.
void updateGradient(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol) override
Update internal risk measure storage for gradient computation.
ROL::Ptr< MixedCVaR< Real > > mqq_
SpectralRisk(ROL::ParameterList &parlist)
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) override
Return risk measure Hessian-times-a-vector.
std::vector< Real > pts_
void resetStorage(UpdateType type) override
Real getValue(const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler) override
Return risk measure value.
SpectralRisk(const std::vector< Real > &pts, const std::vector< Real > &wts, const ROL::Ptr< PlusFunction< Real > > &pf)
void setHessVecStorage(const Ptr< ScalarController< Real > > &gradvec_storage, const Ptr< VectorController< Real > > &hessvec_storage) override
Defines the linear algebra or vector space interface.