ROL
ROL_QuantileQuadrangle.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_SMOOTHCVARQUAD_HPP
11#define ROL_SMOOTHCVARQUAD_HPP
12
14#include "ROL_PlusFunction.hpp"
15
60namespace ROL {
61
62template<class Real>
63class QuantileQuadrangle : public ExpectationQuad<Real> {
64private:
65
66 ROL::Ptr<PlusFunction<Real> > pf_;
67
68 Real prob_;
69 Real lam_;
70 Real eps_;
71
72 Real alpha_;
73 Real beta_;
74
75 void parseParameterList(ROL::ParameterList &parlist) {
76 std::string type = parlist.sublist("SOL").get("Type","Risk Averse");
77 ROL::ParameterList list;
78 if (type == "Risk Averse") {
79 list = parlist.sublist("SOL").sublist("Risk Measure").sublist("CVaR");
80 }
81 else if (type == "Error") {
82 list = parlist.sublist("SOL").sublist("Error Measure").sublist("Koenker-Bassett");
83 }
84 else if (type == "Deviation") {
85 list = parlist.sublist("SOL").sublist("Deviation Measure").sublist("CVaR");
86 }
87 else if (type == "Regret") {
88 list = parlist.sublist("SOL").sublist("Regret Measure").sublist("Mean Absolute Loss");
89 }
90 // Get CVaR parameters
91 prob_ = list.get<Real>("Confidence Level");
92 lam_ = list.get<Real>("Convex Combination Parameter");
93 eps_ = list.get<Real>("Smoothing Parameter");
94 // Build plus function
95 pf_ = ROL::makePtr<PlusFunction<Real>>(list);
96 }
97
98 void checkInputs(void) const {
99 Real zero(0), one(1);
100 ROL_TEST_FOR_EXCEPTION((prob_ <= zero) || (prob_ >= one), std::invalid_argument,
101 ">>> ERROR (ROL::QuantileQuadrangle): Confidence level must be between 0 and 1!");
102 ROL_TEST_FOR_EXCEPTION((lam_ < zero) || (lam_ > one), std::invalid_argument,
103 ">>> ERROR (ROL::QuantileQuadrangle): Convex combination parameter must be positive!");
104 ROL_TEST_FOR_EXCEPTION((eps_ <= zero), std::invalid_argument,
105 ">>> ERROR (ROL::QuantileQuadrangle): Smoothing parameter must be positive!");
106 ROL_TEST_FOR_EXCEPTION(pf_ == ROL::nullPtr, std::invalid_argument,
107 ">>> ERROR (ROL::QuantileQuadrangle): PlusFunction pointer is null!");
108 }
109
110 void setParameters(void) {
111 Real one(1);
112 alpha_ = lam_;
113 beta_ = (one-alpha_*prob_)/(one-prob_);
114 }
115
116public:
123 QuantileQuadrangle(Real prob, Real eps, ROL::Ptr<PlusFunction<Real> > &pf )
124 : ExpectationQuad<Real>(), prob_(prob), lam_(0), eps_(eps), pf_(pf) {
125 checkInputs();
127 }
128
138 QuantileQuadrangle(Real prob, Real lam, Real eps,
139 ROL::Ptr<PlusFunction<Real> > &pf )
140 : ExpectationQuad<Real>(), prob_(prob), lam_(lam), eps_(eps), pf_(pf) {
141 checkInputs();
143 }
144
156 QuantileQuadrangle(ROL::ParameterList &parlist) : ExpectationQuad<Real>() {
157 parseParameterList(parlist);
158 // Check inputs
159 checkInputs();
161 }
162
163 Real error(Real x, int deriv = 0) {
164 Real one(1);
165 Real err = (beta_-one)*pf_->evaluate(x,deriv)
166 + ((deriv%2) ? -one : one)*(one-alpha_)*pf_->evaluate(-x,deriv);
167 return err;
168 }
169
170 Real regret(Real x, int deriv = 0) {
171 Real zero(0), one(1);
172 Real X = ((deriv==0) ? x : ((deriv==1) ? one : zero));
173 Real reg = error(x,deriv) + X;
174 return reg;
175 }
176
177 void check(void) {
179 // Check v'(eps)
180 Real x = eps_, two(2), p1(0.1), one(1), zero(0);
181 Real vx(0), vy(0);
182 Real dv = regret(x,1);
183 Real t(1), diff(0), err(0);
184 std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(eps) is correct? \n";
185 std::cout << std::right << std::setw(20) << "t"
186 << std::setw(20) << "v'(x)"
187 << std::setw(20) << "(v(x+t)-v(x-t))/2t"
188 << std::setw(20) << "Error"
189 << "\n";
190 for (int i = 0; i < 13; i++) {
191 vy = regret(x+t,0);
192 vx = regret(x-t,0);
193 diff = (vy-vx)/(two*t);
194 err = std::abs(diff-dv);
195 std::cout << std::scientific << std::setprecision(11) << std::right
196 << std::setw(20) << t
197 << std::setw(20) << dv
198 << std::setw(20) << diff
199 << std::setw(20) << err
200 << "\n";
201 t *= p1;
202 }
203 std::cout << "\n";
204 // check v''(eps)
205 vx = zero;
206 vy = zero;
207 dv = regret(x,2);
208 t = one;
209 diff = zero;
210 err = zero;
211 std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(eps) is correct? \n";
212 std::cout << std::right << std::setw(20) << "t"
213 << std::setw(20) << "v''(x)"
214 << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
215 << std::setw(20) << "Error"
216 << "\n";
217 for (int i = 0; i < 13; i++) {
218 vy = regret(x+t,1);
219 vx = regret(x-t,1);
220 diff = (vy-vx)/(two*t);
221 err = std::abs(diff-dv);
222 std::cout << std::scientific << std::setprecision(11) << std::right
223 << std::setw(20) << t
224 << std::setw(20) << dv
225 << std::setw(20) << diff
226 << std::setw(20) << err
227 << "\n";
228 t *= p1;
229 }
230 std::cout << "\n";
231 // Check v'(0)
232 x = zero;
233 vx = zero;
234 vy = zero;
235 dv = regret(x,1);
236 t = one;
237 diff = zero;
238 err = zero;
239 std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(0) is correct? \n";
240 std::cout << std::right << std::setw(20) << "t"
241 << std::setw(20) << "v'(x)"
242 << std::setw(20) << "(v(x+t)-v(x-t))/2t"
243 << std::setw(20) << "Error"
244 << "\n";
245 for (int i = 0; i < 13; i++) {
246 vy = regret(x+t,0);
247 vx = regret(x-t,0);
248 diff = (vy-vx)/(two*t);
249 err = std::abs(diff-dv);
250 std::cout << std::scientific << std::setprecision(11) << std::right
251 << std::setw(20) << t
252 << std::setw(20) << dv
253 << std::setw(20) << diff
254 << std::setw(20) << err
255 << "\n";
256 t *= p1;
257 }
258 std::cout << "\n";
259 // check v''(eps)
260 vx = zero;
261 vy = zero;
262 dv = regret(x,2);
263 t = one;
264 diff = zero;
265 err = zero;
266 std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(0) is correct? \n";
267 std::cout << std::right << std::setw(20) << "t"
268 << std::setw(20) << "v''(x)"
269 << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
270 << std::setw(20) << "Error"
271 << "\n";
272 for (int i = 0; i < 13; i++) {
273 vy = regret(x+t,1);
274 vx = regret(x-t,1);
275 diff = (vy-vx)/(two*t);
276 err = std::abs(diff-dv);
277 std::cout << std::scientific << std::setprecision(11) << std::right
278 << std::setw(20) << t
279 << std::setw(20) << dv
280 << std::setw(20) << diff
281 << std::setw(20) << err
282 << "\n";
283 t *= p1;
284 }
285 std::cout << "\n";
286 // Check v'(0)
287 x = -eps_;
288 vx = zero;
289 vy = zero;
290 dv = regret(x,1);
291 t = one;
292 diff = zero;
293 err = zero;
294 std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(-eps) is correct? \n";
295 std::cout << std::right << std::setw(20) << "t"
296 << std::setw(20) << "v'(x)"
297 << std::setw(20) << "(v(x+t)-v(x-t))/2t"
298 << std::setw(20) << "Error"
299 << "\n";
300 for (int i = 0; i < 13; i++) {
301 vy = regret(x+t,0);
302 vx = regret(x-t,0);
303 diff = (vy-vx)/(two*t);
304 err = std::abs(diff-dv);
305 std::cout << std::scientific << std::setprecision(11) << std::right
306 << std::setw(20) << t
307 << std::setw(20) << dv
308 << std::setw(20) << diff
309 << std::setw(20) << err
310 << "\n";
311 t *= p1;
312 }
313 std::cout << "\n";
314 // check v''(eps)
315 vx = zero;
316 vy = zero;
317 dv = regret(x,2);
318 t = one;
319 diff = zero;
320 err = zero;
321 std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(-eps) is correct? \n";
322 std::cout << std::right << std::setw(20) << "t"
323 << std::setw(20) << "v''(x)"
324 << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
325 << std::setw(20) << "Error"
326 << "\n";
327 for (int i = 0; i < 13; i++) {
328 vy = regret(x+t,1);
329 vx = regret(x-t,1);
330 diff = (vy-vx)/(two*t);
331 err = std::abs(diff-dv);
332 std::cout << std::scientific << std::setprecision(11) << std::right
333 << std::setw(20) << t
334 << std::setw(20) << dv
335 << std::setw(20) << diff
336 << std::setw(20) << err
337 << "\n";
338 t *= p1;
339 }
340 std::cout << "\n";
341 }
342
343};
344
345}
346#endif
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Provides a general interface for risk and error measures generated through the expectation risk quadr...
virtual void check(void)
Run default derivative tests for the scalar regret function.
Provides an interface for a convex combination of the expected value and the conditional value-at-ris...
void check(void)
Run default derivative tests for the scalar regret function.
void parseParameterList(ROL::ParameterList &parlist)
Real error(Real x, int deriv=0)
Evaluate the scalar error function at x.
QuantileQuadrangle(Real prob, Real eps, ROL::Ptr< PlusFunction< Real > > &pf)
Constructor.
QuantileQuadrangle(Real prob, Real lam, Real eps, ROL::Ptr< PlusFunction< Real > > &pf)
Constructor.
ROL::Ptr< PlusFunction< Real > > pf_
Real regret(Real x, int deriv=0)
Evaluate the scalar regret function at x.
QuantileQuadrangle(ROL::ParameterList &parlist)
Constructor.