/usr/include/trilinos/ROL_QuantileQuadrangle.hpp is in libtrilinos-rol-dev 12.10.1-3.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 | // @HEADER
// ************************************************************************
//
// Rapid Optimization Library (ROL) Package
// Copyright (2014) Sandia Corporation
//
// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
// license for use of this work by or on behalf of the U.S. Government.
//
// Redistribution and use in source and binary forms, with or without
// modification, are permitted provided that the following conditions are
// met:
//
// 1. Redistributions of source code must retain the above copyright
// notice, this list of conditions and the following disclaimer.
//
// 2. Redistributions in binary form must reproduce the above copyright
// notice, this list of conditions and the following disclaimer in the
// documentation and/or other materials provided with the distribution.
//
// 3. Neither the name of the Corporation nor the names of the
// contributors may be used to endorse or promote products derived from
// this software without specific prior written permission.
//
// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
//
// Questions? Contact lead developers:
// Drew Kouri (dpkouri@sandia.gov) and
// Denis Ridzal (dridzal@sandia.gov)
//
// ************************************************************************
// @HEADER
#ifndef ROL_SMOOTHCVARQUAD_HPP
#define ROL_SMOOTHCVARQUAD_HPP
#include "ROL_ExpectationQuad.hpp"
#include "ROL_PlusFunction.hpp"
/** @ingroup risk_group
\class ROL::QuantileQuadrangle
\brief Provides an interface for a convex combination of the
expected value and the conditional value-at-risk using the
expectation risk quadrangle.
The conditional value-at-risk (also called the average value-at-risk
or the expected shortfall) with confidence level \f$0\le \beta < 1\f$
is
\f[
\mathcal{R}(X) = \inf_{t\in\mathbb{R}} \left\{
t + \frac{1}{1-\beta} \mathbb{E}\left[(X-t)_+\right]
\right\}
\f]
where \f$(x)_+ = \max\{0,x\}\f$. If the distribution of \f$X\f$ is
continuous, then \f$\mathcal{R}\f$ is the conditional expectation of
\f$X\f$ exceeding the \f$\beta\f$-quantile of \f$X\f$ and the optimal
\f$t\f$ is the \f$\beta\f$-quantile.
Additionally, \f$\mathcal{R}\f$ is a law-invariant coherent risk measure.
This class defines a convex combination of expected value and the conditional
value-at-risk using the expectation risk quadrangle. In this case, the scalar
regret function is
\f[
v(x) = \alpha (x)_+ - \lambda (-x)_+
\f]
for \f$\alpha > 1\f$ and \f$0 \le \lambda < 1\f$. The associated confidence
level for the conditional value-at-risk is
\f[
\beta = \frac{\alpha-1}{\alpha-\lambda}.
\f]
This convex combination of expected value and the conditional value-at-risk
is then realized as
\f[
\mathcal{R}(X) = \inf_{t\in\mathbb{R}}\left\{
t + \mathbb{E}[v(X-t)] \right\}.
\f]
ROL implements this by augmenting the optimization vector \f$x_0\f$ with
the parameter \f$t\f$, then minimizes jointly for \f$(x_0,t)\f$.
When using derivative-based optimization, the user can provide a smooth
approximation of \f$(\cdot)_+\f$ using the ROL::PlusFunction class.
*/
namespace ROL {
template<class Real>
class QuantileQuadrangle : public ExpectationQuad<Real> {
private:
Teuchos::RCP<PlusFunction<Real> > pf_;
Real prob_;
Real lam_;
Real eps_;
Real alpha_;
Real beta_;
void checkInputs(void) const {
Real zero(0), one(1);
TEUCHOS_TEST_FOR_EXCEPTION((prob_ <= zero) || (prob_ >= one), std::invalid_argument,
">>> ERROR (ROL::QuantileQuadrangle): Confidence level must be between 0 and 1!");
TEUCHOS_TEST_FOR_EXCEPTION((lam_ < zero) || (lam_ > one), std::invalid_argument,
">>> ERROR (ROL::QuantileQuadrangle): Convex combination parameter must be positive!");
TEUCHOS_TEST_FOR_EXCEPTION((eps_ <= zero), std::invalid_argument,
">>> ERROR (ROL::QuantileQuadrangle): Smoothing parameter must be positive!");
TEUCHOS_TEST_FOR_EXCEPTION(pf_ == Teuchos::null, std::invalid_argument,
">>> ERROR (ROL::QuantileQuadrangle): PlusFunction pointer is null!");
}
void setParameters(void) {
Real one(1);
alpha_ = lam_;
beta_ = (one-alpha_*prob_)/(one-prob_);
}
public:
/** \brief Constructor.
@param[in] prob is the confidence level
@param[in] eps is the smoothing parameter for the plus function approximation
@param[in] pf is the plus function or an approximation
*/
QuantileQuadrangle(Real prob, Real eps, Teuchos::RCP<PlusFunction<Real> > &pf )
: ExpectationQuad<Real>(), prob_(prob), lam_(0), eps_(eps), pf_(pf) {
checkInputs();
setParameters();
}
/** \brief Constructor.
@param[in] prob is the confidence level
@param[in] lam is the convex combination parameter (coeff=0
corresponds to the expected value whereas coeff=1
corresponds to the conditional value-at-risk)
@param[in] eps is the smoothing parameter for the plus function approximation
@param[in] pf is the plus function or an approximation
*/
QuantileQuadrangle(Real prob, Real lam, Real eps,
Teuchos::RCP<PlusFunction<Real> > &pf )
: ExpectationQuad<Real>(), prob_(prob), lam_(lam), eps_(eps), pf_(pf) {
checkInputs();
setParameters();
}
/** \brief Constructor.
@param[in] parlist is a parameter list specifying inputs
parlist should contain sublists "SOL"->"Risk Measure"->"CVaR" and
within the "CVaR" sublist should have the following parameters
\li "Confidence Level" (between 0 and 1)
\li "Convex Combination Parameter" (between 0 and 1)
\li "Smoothing Parameter" (must be positive)
\li A sublist for plus function information.
*/
QuantileQuadrangle(Teuchos::ParameterList &parlist) : ExpectationQuad<Real>() {
Teuchos::ParameterList& list
= parlist.sublist("SOL").sublist("Risk Measure").sublist("Quantile-Based Quadrangle");
// Get CVaR parameters
prob_ = list.get<Real>("Confidence Level");
lam_ = list.get<Real>("Convex Combination Parameter");
eps_ = list.get<Real>("Smoothing Parameter");
// Build plus function
pf_ = Teuchos::rcp( new PlusFunction<Real>(list) );
// Check inputs
checkInputs();
setParameters();
}
Real error(Real x, int deriv = 0) {
Real one(1);
Real err = (beta_-one)*pf_->evaluate(x,deriv)
+ ((deriv%2) ? -one : one)*(one-alpha_)*pf_->evaluate(-x,deriv);
return err;
}
Real regret(Real x, int deriv = 0) {
Real zero(0), one(1);
Real X = ((deriv==0) ? x : ((deriv==1) ? one : zero));
Real reg = error(x,deriv) + X;
return reg;
}
void checkRegret(void) {
ExpectationQuad<Real>::checkRegret();
// Check v'(eps)
Real x = eps_, two(2), p1(0.1), one(1), zero(0);
Real vx(0), vy(0);
Real dv = regret(x,1);
Real t(1), diff(0), err(0);
std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(eps) is correct? \n";
std::cout << std::right << std::setw(20) << "t"
<< std::setw(20) << "v'(x)"
<< std::setw(20) << "(v(x+t)-v(x-t))/2t"
<< std::setw(20) << "Error"
<< "\n";
for (int i = 0; i < 13; i++) {
vy = regret(x+t,0);
vx = regret(x-t,0);
diff = (vy-vx)/(two*t);
err = std::abs(diff-dv);
std::cout << std::scientific << std::setprecision(11) << std::right
<< std::setw(20) << t
<< std::setw(20) << dv
<< std::setw(20) << diff
<< std::setw(20) << err
<< "\n";
t *= p1;
}
std::cout << "\n";
// check v''(eps)
vx = zero;
vy = zero;
dv = regret(x,2);
t = one;
diff = zero;
err = zero;
std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(eps) is correct? \n";
std::cout << std::right << std::setw(20) << "t"
<< std::setw(20) << "v''(x)"
<< std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
<< std::setw(20) << "Error"
<< "\n";
for (int i = 0; i < 13; i++) {
vy = regret(x+t,1);
vx = regret(x-t,1);
diff = (vy-vx)/(two*t);
err = std::abs(diff-dv);
std::cout << std::scientific << std::setprecision(11) << std::right
<< std::setw(20) << t
<< std::setw(20) << dv
<< std::setw(20) << diff
<< std::setw(20) << err
<< "\n";
t *= p1;
}
std::cout << "\n";
// Check v'(0)
x = zero;
vx = zero;
vy = zero;
dv = regret(x,1);
t = one;
diff = zero;
err = zero;
std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(0) is correct? \n";
std::cout << std::right << std::setw(20) << "t"
<< std::setw(20) << "v'(x)"
<< std::setw(20) << "(v(x+t)-v(x-t))/2t"
<< std::setw(20) << "Error"
<< "\n";
for (int i = 0; i < 13; i++) {
vy = regret(x+t,0);
vx = regret(x-t,0);
diff = (vy-vx)/(two*t);
err = std::abs(diff-dv);
std::cout << std::scientific << std::setprecision(11) << std::right
<< std::setw(20) << t
<< std::setw(20) << dv
<< std::setw(20) << diff
<< std::setw(20) << err
<< "\n";
t *= p1;
}
std::cout << "\n";
// check v''(eps)
vx = zero;
vy = zero;
dv = regret(x,2);
t = one;
diff = zero;
err = zero;
std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(0) is correct? \n";
std::cout << std::right << std::setw(20) << "t"
<< std::setw(20) << "v''(x)"
<< std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
<< std::setw(20) << "Error"
<< "\n";
for (int i = 0; i < 13; i++) {
vy = regret(x+t,1);
vx = regret(x-t,1);
diff = (vy-vx)/(two*t);
err = std::abs(diff-dv);
std::cout << std::scientific << std::setprecision(11) << std::right
<< std::setw(20) << t
<< std::setw(20) << dv
<< std::setw(20) << diff
<< std::setw(20) << err
<< "\n";
t *= p1;
}
std::cout << "\n";
// Check v'(0)
x = -eps_;
vx = zero;
vy = zero;
dv = regret(x,1);
t = one;
diff = zero;
err = zero;
std::cout << std::right << std::setw(20) << "CHECK REGRET: v'(-eps) is correct? \n";
std::cout << std::right << std::setw(20) << "t"
<< std::setw(20) << "v'(x)"
<< std::setw(20) << "(v(x+t)-v(x-t))/2t"
<< std::setw(20) << "Error"
<< "\n";
for (int i = 0; i < 13; i++) {
vy = regret(x+t,0);
vx = regret(x-t,0);
diff = (vy-vx)/(two*t);
err = std::abs(diff-dv);
std::cout << std::scientific << std::setprecision(11) << std::right
<< std::setw(20) << t
<< std::setw(20) << dv
<< std::setw(20) << diff
<< std::setw(20) << err
<< "\n";
t *= p1;
}
std::cout << "\n";
// check v''(eps)
vx = zero;
vy = zero;
dv = regret(x,2);
t = one;
diff = zero;
err = zero;
std::cout << std::right << std::setw(20) << "CHECK REGRET: v''(-eps) is correct? \n";
std::cout << std::right << std::setw(20) << "t"
<< std::setw(20) << "v''(x)"
<< std::setw(20) << "(v'(x+t)-v'(x-t))/2t"
<< std::setw(20) << "Error"
<< "\n";
for (int i = 0; i < 13; i++) {
vy = regret(x+t,1);
vx = regret(x-t,1);
diff = (vy-vx)/(two*t);
err = std::abs(diff-dv);
std::cout << std::scientific << std::setprecision(11) << std::right
<< std::setw(20) << t
<< std::setw(20) << dv
<< std::setw(20) << diff
<< std::setw(20) << err
<< "\n";
t *= p1;
}
std::cout << "\n";
}
};
}
#endif
|