ROL
ROL_LogQuantileQuadrangle.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_LOGQUANTILEQUAD_HPP
11#define ROL_LOGQUANTILEQUAD_HPP
12
14#include "ROL_PlusFunction.hpp"
15
39namespace ROL {
40
41template<class Real>
43private:
44
45 ROL::Ptr<PlusFunction<Real> > pf_;
46
47 Real alpha_;
48 Real rate_;
49 Real eps_;
50
51 void parseParameterList(ROL::ParameterList &parlist) {
52 std::string type = parlist.sublist("SOL").get("Type","Risk Averse");
53 ROL::ParameterList list;
54 if (type == "Risk Averse") {
55 list = parlist.sublist("SOL").sublist("Risk Measure").sublist("Log Quantile");
56 }
57 else if (type == "Error") {
58 list = parlist.sublist("SOL").sublist("Error Measure").sublist("Log Quantile");
59 }
60 else if (type == "Deviation") {
61 list = parlist.sublist("SOL").sublist("Deviation Measure").sublist("Log Quantile");
62 }
63 else if (type == "Regret") {
64 list = parlist.sublist("SOL").sublist("Regret Measure").sublist("Log Quantile");
65 }
66 // Check CVaR inputs
67 alpha_ = list.get<Real>("Slope for Linear Growth");
68 rate_ = list.get<Real>("Rate for Exponential Growth");
69 eps_ = list.get<Real>("Smoothing Parameter");
70 // Build plus function
71 pf_ = ROL::makePtr<PlusFunction<Real>>(list);
72 }
73
74 void checkInputs(void) const {
75 Real zero(0), one(1);
76 ROL_TEST_FOR_EXCEPTION((alpha_ < zero) || (alpha_ >= one), std::invalid_argument,
77 ">>> ERROR (ROL::LogQuantileQuadrangle): Linear growth rate must be between 0 and 1!");
78 ROL_TEST_FOR_EXCEPTION((rate_ <= zero), std::invalid_argument,
79 ">>> ERROR (ROL::LogQuantileQuadrangle): Exponential growth rate must be positive!");
80 ROL_TEST_FOR_EXCEPTION((eps_ <= zero), std::invalid_argument,
81 ">>> ERROR (ROL::LogQuantileQuadrangle): Smoothing parameter must be positive!");
82 ROL_TEST_FOR_EXCEPTION(pf_ == ROL::nullPtr, std::invalid_argument,
83 ">>> ERROR (ROL::LogQuantileQuadrangle): PlusFunction pointer is null!");
84 }
85
86public:
94 LogQuantileQuadrangle(Real alpha, Real rate, Real eps,
95 ROL::Ptr<PlusFunction<Real> > &pf )
96 : ExpectationQuad<Real>(), alpha_(alpha), rate_(rate), eps_(eps), pf_(pf) {
98 }
99
111 LogQuantileQuadrangle(ROL::ParameterList &parlist)
112 : ExpectationQuad<Real>() {
113 parseParameterList(parlist);
114 checkInputs();
115 }
116
117 Real error(Real x, int deriv = 0) {
118 Real zero(0), one(1);
119 ROL_TEST_FOR_EXCEPTION( (deriv > 2), std::invalid_argument,
120 ">>> ERROR (ROL::LogQuantileQuadrangle::error): deriv greater than 2!");
121 ROL_TEST_FOR_EXCEPTION( (deriv < 0), std::invalid_argument,
122 ">>> ERROR (ROL::LogQuantileQuadrangle::error): deriv less than 0!");
123
124 Real X = ((deriv == 0) ? x : ((deriv == 1) ? one : zero));
125 return regret(x,deriv) - X;
126 }
127
128 Real regret(Real x, int deriv = 0) {
129 Real zero(0), one(1);
130 ROL_TEST_FOR_EXCEPTION( (deriv > 2), std::invalid_argument,
131 ">>> ERROR (ROL::LogQuantileQuadrangle::regret): deriv greater than 2!");
132 ROL_TEST_FOR_EXCEPTION( (deriv < 0), std::invalid_argument,
133 ">>> ERROR (ROL::LogQuantileQuadrangle::regret): deriv less than 0!");
134
135 Real arg = std::exp(rate_*x);
136 Real sarg = rate_*arg;
137 Real reg = (pf_->evaluate(arg-one,deriv) *
138 ((deriv == 0) ? one/rate_ : ((deriv == 1) ? arg : sarg*arg))
139 + ((deriv == 2) ? pf_->evaluate(arg-one,deriv-1)*sarg : zero))
140 + ((deriv%2 == 0) ? -one : one) * alpha_ * pf_->evaluate(-x,deriv);
141 return reg;
142 }
143
144 void check(void) {
146 // Check v'(eps)
147 Real x = eps_, two(2), p1(0.1), zero(0), one(1);
148 Real vx(0), vy(0);
149 Real dv = regret(x,1);
150 Real t(1), diff(0), err(0);
151 std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(eps) is correct? \n";
152 std::cout << std::right << std::setw(20) << "t"
153 << std::setw(20) << "v'(x)"
154 << std::setw(20) << "(v(x+t)-v(x-t))/2t"
155 << std::setw(20) << "Error"
156 << "\n";
157 for (int i = 0; i < 13; i++) {
158 vy = regret(x+t,0);
159 vx = regret(x-t,0);
160 diff = (vy-vx)/(two*t);
161 err = std::abs(diff-dv);
162 std::cout << std::scientific << std::setprecision(11) << std::right
163 << std::setw(20) << t
164 << std::setw(20) << dv
165 << std::setw(20) << diff
166 << std::setw(20) << err
167 << "\n";
168 t *= p1;
169 }
170 std::cout << "\n";
171 // check v''(eps)
172 vx = zero;
173 vy = zero;
174 dv = regret(x,2);
175 t = one;
176 diff = zero;
177 err = zero;
178 std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(eps) is correct? \n";
179 std::cout << std::right << std::setw(20) << "t"
180 << std::setw(20) << "v''(x)"
181 << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
182 << std::setw(20) << "Error"
183 << "\n";
184 for (int i = 0; i < 13; i++) {
185 vy = regret(x+t,1);
186 vx = regret(x-t,1);
187 diff = (vy-vx)/(two*t);
188 err = std::abs(diff-dv);
189 std::cout << std::scientific << std::setprecision(11) << std::right
190 << std::setw(20) << t
191 << std::setw(20) << dv
192 << std::setw(20) << diff
193 << std::setw(20) << err
194 << "\n";
195 t *= p1;
196 }
197 std::cout << "\n";
198 // Check v'(0)
199 x = zero;
200 vx = zero;
201 vy = zero;
202 dv = regret(x,1);
203 t = one;
204 diff = zero;
205 err = zero;
206 std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(0) is correct? \n";
207 std::cout << std::right << std::setw(20) << "t"
208 << std::setw(20) << "v'(x)"
209 << std::setw(20) << "(v(x+t)-v(x-t))/2t"
210 << std::setw(20) << "Error"
211 << "\n";
212 for (int i = 0; i < 13; i++) {
213 vy = regret(x+t,0);
214 vx = regret(x-t,0);
215 diff = (vy-vx)/(two*t);
216 err = std::abs(diff-dv);
217 std::cout << std::scientific << std::setprecision(11) << std::right
218 << std::setw(20) << t
219 << std::setw(20) << dv
220 << std::setw(20) << diff
221 << std::setw(20) << err
222 << "\n";
223 t *= p1;
224 }
225 std::cout << "\n";
226 // check v''(eps)
227 vx = zero;
228 vy = zero;
229 dv = regret(x,2);
230 t = one;
231 diff = zero;
232 err = zero;
233 std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(0) is correct? \n";
234 std::cout << std::right << std::setw(20) << "t"
235 << std::setw(20) << "v''(x)"
236 << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
237 << std::setw(20) << "Error"
238 << "\n";
239 for (int i = 0; i < 13; i++) {
240 vy = regret(x+t,1);
241 vx = regret(x-t,1);
242 diff = (vy-vx)/(two*t);
243 err = std::abs(diff-dv);
244 std::cout << std::scientific << std::setprecision(11) << std::right
245 << std::setw(20) << t
246 << std::setw(20) << dv
247 << std::setw(20) << diff
248 << std::setw(20) << err
249 << "\n";
250 t *= p1;
251 }
252 std::cout << "\n";
253 // Check v'(0)
254 x = -eps_;
255 vx = zero;
256 vy = zero;
257 dv = regret(x,1);
258 t = one;
259 diff = zero;
260 err = zero;
261 std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(-eps) is correct? \n";
262 std::cout << std::right << std::setw(20) << "t"
263 << std::setw(20) << "v'(x)"
264 << std::setw(20) << "(v(x+t)-v(x-t))/2t"
265 << std::setw(20) << "Error"
266 << "\n";
267 for (int i = 0; i < 13; i++) {
268 vy = regret(x+t,0);
269 vx = regret(x-t,0);
270 diff = (vy-vx)/(two*t);
271 err = std::abs(diff-dv);
272 std::cout << std::scientific << std::setprecision(11) << std::right
273 << std::setw(20) << t
274 << std::setw(20) << dv
275 << std::setw(20) << diff
276 << std::setw(20) << err
277 << "\n";
278 t *= p1;
279 }
280 std::cout << "\n";
281 // check v''(eps)
282 vx = zero;
283 vy = zero;
284 dv = regret(x,2);
285 t = one;
286 diff = zero;
287 err = zero;
288 std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(-eps) is correct? \n";
289 std::cout << std::right << std::setw(20) << "t"
290 << std::setw(20) << "v''(x)"
291 << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
292 << std::setw(20) << "Error"
293 << "\n";
294 for (int i = 0; i < 13; i++) {
295 vy = regret(x+t,1);
296 vx = regret(x-t,1);
297 diff = (vy-vx)/(two*t);
298 err = std::abs(diff-dv);
299 std::cout << std::scientific << std::setprecision(11) << std::right
300 << std::setw(20) << t
301 << std::setw(20) << dv
302 << std::setw(20) << diff
303 << std::setw(20) << err
304 << "\n";
305 t *= p1;
306 }
307 std::cout << "\n";
308 }
309
310};
311
312}
313#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 the conditioanl entropic risk using the expectation risk quadrangle.
void check(void)
Run default derivative tests for the scalar regret function.
Real error(Real x, int deriv=0)
Evaluate the scalar error function at x.
void parseParameterList(ROL::ParameterList &parlist)
LogQuantileQuadrangle(Real alpha, Real rate, Real eps, ROL::Ptr< PlusFunction< Real > > &pf)
Constructor.
LogQuantileQuadrangle(ROL::ParameterList &parlist)
Constructor.
ROL::Ptr< PlusFunction< Real > > pf_
Real regret(Real x, int deriv=0)
Evaluate the scalar regret function at x.