ROL
ROL_QuantileQuadrangle.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_SMOOTHCVARQUAD_HPP
45 #define ROL_SMOOTHCVARQUAD_HPP
46 
47 #include "ROL_ExpectationQuad.hpp"
48 #include "ROL_PlusFunction.hpp"
49 
94 namespace ROL {
95 
96 template<class Real>
97 class QuantileQuadrangle : public ExpectationQuad<Real> {
98 private:
99 
100  ROL::Ptr<PlusFunction<Real> > pf_;
101 
102  Real prob_;
103  Real lam_;
104  Real eps_;
105 
106  Real alpha_;
107  Real beta_;
108 
109  void parseParameterList(ROL::ParameterList &parlist) {
110  std::string type = parlist.sublist("SOL").get("Stochastic Component Type","Risk Averse");
111  ROL::ParameterList list;
112  if (type == "Risk Averse") {
113  list = parlist.sublist("SOL").sublist("Risk Measure").sublist("CVaR");
114  }
115  else if (type == "Error") {
116  list = parlist.sublist("SOL").sublist("Error Measure").sublist("Koenker-Bassett");
117  }
118  else if (type == "Deviation") {
119  list = parlist.sublist("SOL").sublist("Deviation Measure").sublist("CVaR");
120  }
121  else if (type == "Regret") {
122  list = parlist.sublist("SOL").sublist("Regret Measure").sublist("Mean Absolute Loss");
123  }
124  // Get CVaR parameters
125  prob_ = list.get<Real>("Confidence Level");
126  lam_ = list.get<Real>("Convex Combination Parameter");
127  eps_ = list.get<Real>("Smoothing Parameter");
128  // Build plus function
129  pf_ = ROL::makePtr<PlusFunction<Real>>(list);
130  }
131 
132  void checkInputs(void) const {
133  Real zero(0), one(1);
134  ROL_TEST_FOR_EXCEPTION((prob_ <= zero) || (prob_ >= one), std::invalid_argument,
135  ">>> ERROR (ROL::QuantileQuadrangle): Confidence level must be between 0 and 1!");
136  ROL_TEST_FOR_EXCEPTION((lam_ < zero) || (lam_ > one), std::invalid_argument,
137  ">>> ERROR (ROL::QuantileQuadrangle): Convex combination parameter must be positive!");
138  ROL_TEST_FOR_EXCEPTION((eps_ <= zero), std::invalid_argument,
139  ">>> ERROR (ROL::QuantileQuadrangle): Smoothing parameter must be positive!");
140  ROL_TEST_FOR_EXCEPTION(pf_ == ROL::nullPtr, std::invalid_argument,
141  ">>> ERROR (ROL::QuantileQuadrangle): PlusFunction pointer is null!");
142  }
143 
144  void setParameters(void) {
145  Real one(1);
146  alpha_ = lam_;
147  beta_ = (one-alpha_*prob_)/(one-prob_);
148  }
149 
150 public:
157  QuantileQuadrangle(Real prob, Real eps, ROL::Ptr<PlusFunction<Real> > &pf )
158  : ExpectationQuad<Real>(), prob_(prob), lam_(0), eps_(eps), pf_(pf) {
159  checkInputs();
160  setParameters();
161  }
162 
172  QuantileQuadrangle(Real prob, Real lam, Real eps,
173  ROL::Ptr<PlusFunction<Real> > &pf )
174  : ExpectationQuad<Real>(), prob_(prob), lam_(lam), eps_(eps), pf_(pf) {
175  checkInputs();
176  setParameters();
177  }
178 
190  QuantileQuadrangle(ROL::ParameterList &parlist) : ExpectationQuad<Real>() {
191  parseParameterList(parlist);
192  // Check inputs
193  checkInputs();
194  setParameters();
195  }
196 
197  Real error(Real x, int deriv = 0) {
198  Real one(1);
199  Real err = (beta_-one)*pf_->evaluate(x,deriv)
200  + ((deriv%2) ? -one : one)*(one-alpha_)*pf_->evaluate(-x,deriv);
201  return err;
202  }
203 
204  Real regret(Real x, int deriv = 0) {
205  Real zero(0), one(1);
206  Real X = ((deriv==0) ? x : ((deriv==1) ? one : zero));
207  Real reg = error(x,deriv) + X;
208  return reg;
209  }
210 
211  void check(void) {
213  // Check v'(eps)
214  Real x = eps_, two(2), p1(0.1), one(1), zero(0);
215  Real vx(0), vy(0);
216  Real dv = regret(x,1);
217  Real t(1), diff(0), err(0);
218  std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(eps) is correct? \n";
219  std::cout << std::right << std::setw(20) << "t"
220  << std::setw(20) << "v'(x)"
221  << std::setw(20) << "(v(x+t)-v(x-t))/2t"
222  << std::setw(20) << "Error"
223  << "\n";
224  for (int i = 0; i < 13; i++) {
225  vy = regret(x+t,0);
226  vx = regret(x-t,0);
227  diff = (vy-vx)/(two*t);
228  err = std::abs(diff-dv);
229  std::cout << std::scientific << std::setprecision(11) << std::right
230  << std::setw(20) << t
231  << std::setw(20) << dv
232  << std::setw(20) << diff
233  << std::setw(20) << err
234  << "\n";
235  t *= p1;
236  }
237  std::cout << "\n";
238  // check v''(eps)
239  vx = zero;
240  vy = zero;
241  dv = regret(x,2);
242  t = one;
243  diff = zero;
244  err = zero;
245  std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(eps) is correct? \n";
246  std::cout << std::right << std::setw(20) << "t"
247  << std::setw(20) << "v''(x)"
248  << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
249  << std::setw(20) << "Error"
250  << "\n";
251  for (int i = 0; i < 13; i++) {
252  vy = regret(x+t,1);
253  vx = regret(x-t,1);
254  diff = (vy-vx)/(two*t);
255  err = std::abs(diff-dv);
256  std::cout << std::scientific << std::setprecision(11) << std::right
257  << std::setw(20) << t
258  << std::setw(20) << dv
259  << std::setw(20) << diff
260  << std::setw(20) << err
261  << "\n";
262  t *= p1;
263  }
264  std::cout << "\n";
265  // Check v'(0)
266  x = zero;
267  vx = zero;
268  vy = zero;
269  dv = regret(x,1);
270  t = one;
271  diff = zero;
272  err = zero;
273  std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(0) is correct? \n";
274  std::cout << std::right << std::setw(20) << "t"
275  << std::setw(20) << "v'(x)"
276  << std::setw(20) << "(v(x+t)-v(x-t))/2t"
277  << std::setw(20) << "Error"
278  << "\n";
279  for (int i = 0; i < 13; i++) {
280  vy = regret(x+t,0);
281  vx = regret(x-t,0);
282  diff = (vy-vx)/(two*t);
283  err = std::abs(diff-dv);
284  std::cout << std::scientific << std::setprecision(11) << std::right
285  << std::setw(20) << t
286  << std::setw(20) << dv
287  << std::setw(20) << diff
288  << std::setw(20) << err
289  << "\n";
290  t *= p1;
291  }
292  std::cout << "\n";
293  // check v''(eps)
294  vx = zero;
295  vy = zero;
296  dv = regret(x,2);
297  t = one;
298  diff = zero;
299  err = zero;
300  std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(0) is correct? \n";
301  std::cout << std::right << std::setw(20) << "t"
302  << std::setw(20) << "v''(x)"
303  << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
304  << std::setw(20) << "Error"
305  << "\n";
306  for (int i = 0; i < 13; i++) {
307  vy = regret(x+t,1);
308  vx = regret(x-t,1);
309  diff = (vy-vx)/(two*t);
310  err = std::abs(diff-dv);
311  std::cout << std::scientific << std::setprecision(11) << std::right
312  << std::setw(20) << t
313  << std::setw(20) << dv
314  << std::setw(20) << diff
315  << std::setw(20) << err
316  << "\n";
317  t *= p1;
318  }
319  std::cout << "\n";
320  // Check v'(0)
321  x = -eps_;
322  vx = zero;
323  vy = zero;
324  dv = regret(x,1);
325  t = one;
326  diff = zero;
327  err = zero;
328  std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(-eps) is correct? \n";
329  std::cout << std::right << std::setw(20) << "t"
330  << std::setw(20) << "v'(x)"
331  << std::setw(20) << "(v(x+t)-v(x-t))/2t"
332  << std::setw(20) << "Error"
333  << "\n";
334  for (int i = 0; i < 13; i++) {
335  vy = regret(x+t,0);
336  vx = regret(x-t,0);
337  diff = (vy-vx)/(two*t);
338  err = std::abs(diff-dv);
339  std::cout << std::scientific << std::setprecision(11) << std::right
340  << std::setw(20) << t
341  << std::setw(20) << dv
342  << std::setw(20) << diff
343  << std::setw(20) << err
344  << "\n";
345  t *= p1;
346  }
347  std::cout << "\n";
348  // check v''(eps)
349  vx = zero;
350  vy = zero;
351  dv = regret(x,2);
352  t = one;
353  diff = zero;
354  err = zero;
355  std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(-eps) is correct? \n";
356  std::cout << std::right << std::setw(20) << "t"
357  << std::setw(20) << "v''(x)"
358  << std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
359  << std::setw(20) << "Error"
360  << "\n";
361  for (int i = 0; i < 13; i++) {
362  vy = regret(x+t,1);
363  vx = regret(x-t,1);
364  diff = (vy-vx)/(two*t);
365  err = std::abs(diff-dv);
366  std::cout << std::scientific << std::setprecision(11) << std::right
367  << std::setw(20) << t
368  << std::setw(20) << dv
369  << std::setw(20) << diff
370  << std::setw(20) << err
371  << "\n";
372  t *= p1;
373  }
374  std::cout << "\n";
375  }
376 
377 };
378 
379 }
380 #endif
QuantileQuadrangle(Real prob, Real eps, ROL::Ptr< PlusFunction< Real > > &pf)
Constructor.
Provides a general interface for risk and error measures generated through the expectation risk quadr...
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...
virtual void check(void)
Run default derivative tests for the scalar regret function.
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
QuantileQuadrangle(ROL::ParameterList &parlist)
Constructor.
Real regret(Real x, int deriv=0)
Evaluate the scalar regret function at x.
QuantileQuadrangle(Real prob, Real lam, Real eps, ROL::Ptr< PlusFunction< Real > > &pf)
Constructor.
Real error(Real x, int deriv=0)
Evaluate the scalar error function at x.
void parseParameterList(ROL::ParameterList &parlist)
ROL::Ptr< PlusFunction< Real > > pf_