ROL
ROL_LogQuantileQuadrangle.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_LOGQUANTILEQUAD_HPP
45 #define ROL_LOGQUANTILEQUAD_HPP
46 
47 #include "ROL_ExpectationQuad.hpp"
48 #include "ROL_PlusFunction.hpp"
49 
73 namespace ROL {
74 
75 template<class Real>
76 class LogQuantileQuadrangle : public ExpectationQuad<Real> {
77 private:
78 
79  ROL::Ptr<PlusFunction<Real> > pf_;
80 
81  Real alpha_;
82  Real rate_;
83  Real eps_;
84 
85  void parseParameterList(ROL::ParameterList &parlist) {
86  std::string type = parlist.sublist("SOL").get("Stochastic Component Type","Risk Averse");
87  ROL::ParameterList list;
88  if (type == "Risk Averse") {
89  list = parlist.sublist("SOL").sublist("Risk Measure").sublist("Log Quantile");
90  }
91  else if (type == "Error") {
92  list = parlist.sublist("SOL").sublist("Error Measure").sublist("Log Quantile");
93  }
94  else if (type == "Deviation") {
95  list = parlist.sublist("SOL").sublist("Deviation Measure").sublist("Log Quantile");
96  }
97  else if (type == "Regret") {
98  list = parlist.sublist("SOL").sublist("Regret Measure").sublist("Log Quantile");
99  }
100  // Check CVaR inputs
101  alpha_ = list.get<Real>("Slope for Linear Growth");
102  rate_ = list.get<Real>("Rate for Exponential Growth");
103  eps_ = list.get<Real>("Smoothing Parameter");
104  // Build plus function
105  pf_ = ROL::makePtr<PlusFunction<Real>>(list);
106  }
107 
108  void checkInputs(void) const {
109  Real zero(0), one(1);
110  ROL_TEST_FOR_EXCEPTION((alpha_ < zero) || (alpha_ >= one), std::invalid_argument,
111  ">>> ERROR (ROL::LogQuantileQuadrangle): Linear growth rate must be between 0 and 1!");
112  ROL_TEST_FOR_EXCEPTION((rate_ <= zero), std::invalid_argument,
113  ">>> ERROR (ROL::LogQuantileQuadrangle): Exponential growth rate must be positive!");
114  ROL_TEST_FOR_EXCEPTION((eps_ <= zero), std::invalid_argument,
115  ">>> ERROR (ROL::LogQuantileQuadrangle): Smoothing parameter must be positive!");
116  ROL_TEST_FOR_EXCEPTION(pf_ == ROL::nullPtr, std::invalid_argument,
117  ">>> ERROR (ROL::LogQuantileQuadrangle): PlusFunction pointer is null!");
118  }
119 
120 public:
128  LogQuantileQuadrangle(Real alpha, Real rate, Real eps,
129  ROL::Ptr<PlusFunction<Real> > &pf )
130  : ExpectationQuad<Real>(), alpha_(alpha), rate_(rate), eps_(eps), pf_(pf) {
131  checkInputs();
132  }
133 
145  LogQuantileQuadrangle(ROL::ParameterList &parlist)
146  : ExpectationQuad<Real>() {
147  parseParameterList(parlist);
148  checkInputs();
149  }
150 
151  Real error(Real x, int deriv = 0) {
152  Real zero(0), one(1);
153  ROL_TEST_FOR_EXCEPTION( (deriv > 2), std::invalid_argument,
154  ">>> ERROR (ROL::LogQuantileQuadrangle::error): deriv greater than 2!");
155  ROL_TEST_FOR_EXCEPTION( (deriv < 0), std::invalid_argument,
156  ">>> ERROR (ROL::LogQuantileQuadrangle::error): deriv less than 0!");
157 
158  Real X = ((deriv == 0) ? x : ((deriv == 1) ? one : zero));
159  return regret(x,deriv) - X;
160  }
161 
162  Real regret(Real x, int deriv = 0) {
163  Real zero(0), one(1);
164  ROL_TEST_FOR_EXCEPTION( (deriv > 2), std::invalid_argument,
165  ">>> ERROR (ROL::LogQuantileQuadrangle::regret): deriv greater than 2!");
166  ROL_TEST_FOR_EXCEPTION( (deriv < 0), std::invalid_argument,
167  ">>> ERROR (ROL::LogQuantileQuadrangle::regret): deriv less than 0!");
168 
169  Real arg = std::exp(rate_*x);
170  Real sarg = rate_*arg;
171  Real reg = (pf_->evaluate(arg-one,deriv) *
172  ((deriv == 0) ? one/rate_ : ((deriv == 1) ? arg : sarg*arg))
173  + ((deriv == 2) ? pf_->evaluate(arg-one,deriv-1)*sarg : zero))
174  + ((deriv%2 == 0) ? -one : one) * alpha_ * pf_->evaluate(-x,deriv);
175  return reg;
176  }
177 
178  void check(void) {
180  // Check v'(eps)
181  Real x = eps_, two(2), p1(0.1), zero(0), one(1);
182  Real vx(0), vy(0);
183  Real dv = regret(x,1);
184  Real t(1), diff(0), err(0);
185  std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(eps) is correct? \n";
186  std::cout << std::right << std::setw(20) << "t"
187  << std::setw(20) << "v'(x)"
188  << std::setw(20) << "(v(x+t)-v(x-t))/2t"
189  << std::setw(20) << "Error"
190  << "\n";
191  for (int i = 0; i < 13; i++) {
192  vy = regret(x+t,0);
193  vx = regret(x-t,0);
194  diff = (vy-vx)/(two*t);
195  err = std::abs(diff-dv);
196  std::cout << std::scientific << std::setprecision(11) << std::right
197  << std::setw(20) << t
198  << std::setw(20) << dv
199  << std::setw(20) << diff
200  << std::setw(20) << err
201  << "\n";
202  t *= p1;
203  }
204  std::cout << "\n";
205  // check v''(eps)
206  vx = zero;
207  vy = zero;
208  dv = regret(x,2);
209  t = one;
210  diff = zero;
211  err = zero;
212  std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(eps) is correct? \n";
213  std::cout << std::right << std::setw(20) << "t"
214  << std::setw(20) << "v''(x)"
215  << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
216  << std::setw(20) << "Error"
217  << "\n";
218  for (int i = 0; i < 13; i++) {
219  vy = regret(x+t,1);
220  vx = regret(x-t,1);
221  diff = (vy-vx)/(two*t);
222  err = std::abs(diff-dv);
223  std::cout << std::scientific << std::setprecision(11) << std::right
224  << std::setw(20) << t
225  << std::setw(20) << dv
226  << std::setw(20) << diff
227  << std::setw(20) << err
228  << "\n";
229  t *= p1;
230  }
231  std::cout << "\n";
232  // Check v'(0)
233  x = zero;
234  vx = zero;
235  vy = zero;
236  dv = regret(x,1);
237  t = one;
238  diff = zero;
239  err = zero;
240  std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(0) is correct? \n";
241  std::cout << std::right << std::setw(20) << "t"
242  << std::setw(20) << "v'(x)"
243  << std::setw(20) << "(v(x+t)-v(x-t))/2t"
244  << std::setw(20) << "Error"
245  << "\n";
246  for (int i = 0; i < 13; i++) {
247  vy = regret(x+t,0);
248  vx = regret(x-t,0);
249  diff = (vy-vx)/(two*t);
250  err = std::abs(diff-dv);
251  std::cout << std::scientific << std::setprecision(11) << std::right
252  << std::setw(20) << t
253  << std::setw(20) << dv
254  << std::setw(20) << diff
255  << std::setw(20) << err
256  << "\n";
257  t *= p1;
258  }
259  std::cout << "\n";
260  // check v''(eps)
261  vx = zero;
262  vy = zero;
263  dv = regret(x,2);
264  t = one;
265  diff = zero;
266  err = zero;
267  std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(0) is correct? \n";
268  std::cout << std::right << std::setw(20) << "t"
269  << std::setw(20) << "v''(x)"
270  << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
271  << std::setw(20) << "Error"
272  << "\n";
273  for (int i = 0; i < 13; i++) {
274  vy = regret(x+t,1);
275  vx = regret(x-t,1);
276  diff = (vy-vx)/(two*t);
277  err = std::abs(diff-dv);
278  std::cout << std::scientific << std::setprecision(11) << std::right
279  << std::setw(20) << t
280  << std::setw(20) << dv
281  << std::setw(20) << diff
282  << std::setw(20) << err
283  << "\n";
284  t *= p1;
285  }
286  std::cout << "\n";
287  // Check v'(0)
288  x = -eps_;
289  vx = zero;
290  vy = zero;
291  dv = regret(x,1);
292  t = one;
293  diff = zero;
294  err = zero;
295  std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(-eps) is correct? \n";
296  std::cout << std::right << std::setw(20) << "t"
297  << std::setw(20) << "v'(x)"
298  << std::setw(20) << "(v(x+t)-v(x-t))/2t"
299  << std::setw(20) << "Error"
300  << "\n";
301  for (int i = 0; i < 13; i++) {
302  vy = regret(x+t,0);
303  vx = regret(x-t,0);
304  diff = (vy-vx)/(two*t);
305  err = std::abs(diff-dv);
306  std::cout << std::scientific << std::setprecision(11) << std::right
307  << std::setw(20) << t
308  << std::setw(20) << dv
309  << std::setw(20) << diff
310  << std::setw(20) << err
311  << "\n";
312  t *= p1;
313  }
314  std::cout << "\n";
315  // check v''(eps)
316  vx = zero;
317  vy = zero;
318  dv = regret(x,2);
319  t = one;
320  diff = zero;
321  err = zero;
322  std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(-eps) is correct? \n";
323  std::cout << std::right << std::setw(20) << "t"
324  << std::setw(20) << "v''(x)"
325  << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
326  << std::setw(20) << "Error"
327  << "\n";
328  for (int i = 0; i < 13; i++) {
329  vy = regret(x+t,1);
330  vx = regret(x-t,1);
331  diff = (vy-vx)/(two*t);
332  err = std::abs(diff-dv);
333  std::cout << std::scientific << std::setprecision(11) << std::right
334  << std::setw(20) << t
335  << std::setw(20) << dv
336  << std::setw(20) << diff
337  << std::setw(20) << err
338  << "\n";
339  t *= p1;
340  }
341  std::cout << "\n";
342  }
343 
344 };
345 
346 }
347 #endif
Real error(Real x, int deriv=0)
Evaluate the scalar error function at x.
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.
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Provides an interface for the conditioanl entropic risk using the expectation risk quadrangle...
LogQuantileQuadrangle(ROL::ParameterList &parlist)
Constructor.
void parseParameterList(ROL::ParameterList &parlist)
ROL::Ptr< PlusFunction< Real > > pf_
void check(void)
Run default derivative tests for the scalar regret function.
Real regret(Real x, int deriv=0)
Evaluate the scalar regret function at x.
LogQuantileQuadrangle(Real alpha, Real rate, Real eps, ROL::Ptr< PlusFunction< Real > > &pf)
Constructor.