This file is indexed.

/usr/include/trilinos/ROL_SpectralRisk.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
// @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_SPECTRALRISK_HPP
#define ROL_SPECTRALRISK_HPP

#include "ROL_MixedQuantileQuadrangle.hpp"
#include "ROL_DistributionFactory.hpp"

/** @ingroup risk_group
    \class ROL::SpectralRisk
    \brief Provides an interface for spectral risk measures.

    Kusuoka's representation for law-invariant risk measures is
    \f[
       \mathcal{R}(X) = \sup_{\mu\in\mathfrak{M}}
          \int_0^1 \mathrm{CVaR}_{\alpha}(X)\,\mathrm{d}\mu(\alpha)
    \f]
    where the conditional value-at-risk (CVaR) with confidence level
    \f$0\le \alpha < 1\f$ is
    \f[
       \mathrm{CVaR}_\alpha(X) = \inf_{t\in\mathbb{R}} \left\{
         t + \frac{1}{1-\alpha} \mathbb{E}\left[(X-t)_+\right]
         \right\}, \quad (x)_+ = \max\{0,x\},
    \f]
    and \f$\mathfrak{M}\f$ is a subset of distributions on the interval
    \f$[0,1)\f$.  By spectral risk measures, we refer to the case where the set
    \f$\mathfrak{M}\f$ is a singleton.  If the distribution
    \f$\mu\in\mathfrak{M}\f$ is discrete, then the corresponding risk measure
    is a mixed quantile quadrangle risk measure.

    If the distribution of \f$X\f$ is continuous, then
    \f$\mathrm{CVaR}_{\alpha}(X)\f$ is the conditional
    expectation of \f$X\f$ exceeding the \f$\alpha\f$-quantile of \f$X\f$ and
    the optimal \f$t\f$ is the \f$\alpha\f$-quantile.
    Additionally, \f$\mathcal{R}\f$ is a law-invariant coherent risk measure.

    ROL implements \f$\mathcal{R}\f$ by approximating the integral with
    Gauss-Chebyshev quadrature of the first kind.  The corresponding quadrature
    points and weights are then used to construct a
    ROL::MixedQuantileQuadrangle risk measure.
    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 SpectralRisk : public RiskMeasure<Real> {
private:
  Teuchos::RCP<MixedQuantileQuadrangle<Real> > mqq_;
  Teuchos::RCP<PlusFunction<Real> > plusFunction_;

  std::vector<Real> wts_;
  std::vector<Real> pts_;

  void checkInputs(Teuchos::RCP<Distribution<Real> > &dist = Teuchos::null) const {
    TEUCHOS_TEST_FOR_EXCEPTION(plusFunction_ == Teuchos::null, std::invalid_argument,
      ">>> ERROR (ROL::SpectralRisk): PlusFunction pointer is null!");
    if ( dist != Teuchos::null) {
      Real lb = dist->lowerBound();
      Real ub = dist->upperBound();
      TEUCHOS_TEST_FOR_EXCEPTION(lb < static_cast<Real>(0), std::invalid_argument,
        ">>> ERROR (ROL::SpectralRisk): Distribution lower bound less than zero!");
      TEUCHOS_TEST_FOR_EXCEPTION(ub > static_cast<Real>(1), std::invalid_argument,
        ">>> ERROR (ROL::SpectralRisk): Distribution upper bound greater than one!");
    }
  }

protected:
  void buildMixedQuantile(const std::vector<Real> &pts, const std::vector<Real> &wts,
                          const Teuchos::RCP<PlusFunction<Real> > &pf) {
     pts_.clear(); pts_.assign(pts.begin(),pts.end());
     wts_.clear(); wts_.assign(wts.begin(),wts.end());
     plusFunction_ = pf;
     mqq_ = Teuchos::rcp(new MixedQuantileQuadrangle<Real>(pts,wts,pf));
  }

  void buildQuadFromDist(std::vector<Real> &pts, std::vector<Real> &wts,
                   const int nQuad, const Teuchos::RCP<Distribution<Real> > &dist) const {
    const Real lo = dist->lowerBound(), hi = dist->upperBound();
    const Real half(0.5), one(1), N(nQuad);
    wts.clear(); wts.resize(nQuad);
    pts.clear(); pts.resize(nQuad);
    if ( hi >= one ) {
      wts[0] = half/(N-half);
      pts[0] = lo;
      for (int i = 1; i < nQuad; ++i) {
        wts[i] = one/(N-half);
        pts[i] = dist->invertCDF(static_cast<Real>(i)/N);
      }
    }
    else {
      wts[0] = half/(N-one);
      pts[0] = lo;
      for (int i = 1; i < nQuad-1; ++i) {
        wts[i] = one/(N-one);
        pts[i] = dist->invertCDF(static_cast<Real>(i)/N);
      }
      wts[nQuad-1] = half/(N-one);
      pts[nQuad-1] = hi;
    }
  }

  void printQuad(const std::vector<Real> &pts,
                 const std::vector<Real> &wts,
                 const bool print = false) const {
    if ( print ) {
      const int nQuad = wts.size();
      std::cout << std::endl;
      std::cout << std::scientific << std::setprecision(15);
      std::cout << std::setw(25) << std::left << "Points"
                << std::setw(25) << std::left << "Weights"
                << std::endl;
      for (int i = 0; i < nQuad; ++i) {
        std::cout << std::setw(25) << std::left << pts[i]
                  << std::setw(25) << std::left << wts[i]
                  << std::endl;
      }
      std::cout << std::endl;
    }
  }

public:
  SpectralRisk(void) : RiskMeasure<Real>() {}

  SpectralRisk( const Teuchos::RCP<Distribution<Real> > &dist,
                const int nQuad,
                const Teuchos::RCP<PlusFunction<Real> > &pf)
    : RiskMeasure<Real>() {
    // Build generalized trapezoidal rule
    std::vector<Real> wts(nQuad), pts(nQuad);
    buildQuadFromDist(pts,wts,nQuad,dist);
    // Build mixed quantile quadrangle risk measure
    buildMixedQuantile(pts,wts,pf);
    // Check inputs
    checkInputs(dist);
  }

  SpectralRisk(Teuchos::ParameterList &parlist)
    : RiskMeasure<Real>() {
    // Parse parameter list
    Teuchos::ParameterList &list
      = parlist.sublist("SOL").sublist("Risk Measure").sublist("Spectral Risk");
    int nQuad  = list.get("Number of Quadrature Points",5);
    bool print = list.get("Print Quadrature to Screen",false);
    // Build distribution
    Teuchos::RCP<Distribution<Real> > dist = DistributionFactory<Real>(list);
    // Build plus function approximation
    Teuchos::RCP<PlusFunction<Real> > pf = Teuchos::rcp(new PlusFunction<Real>(list));
    // Build generalized trapezoidal rule
    std::vector<Real> wts(nQuad), pts(nQuad);
    buildQuadFromDist(pts,wts,nQuad,dist);
    printQuad(pts,wts,print);
    // Build mixed quantile quadrangle risk measure
    buildMixedQuantile(pts,wts,pf);
    // Check inputs
    checkInputs(dist);
  }

  SpectralRisk( const std::vector<Real> &pts, const std::vector<Real> &wts,
                const Teuchos::RCP<PlusFunction<Real> > &pf)
    : RiskMeasure<Real>() {
    buildMixedQuantile(pts,wts,pf);
    // Check inputs
    checkInputs();
  }

  Real computeStatistic(const Vector<Real> &x) const {
    std::vector<Real> xstat;
    Teuchos::dyn_cast<const RiskVector<Real> >(x).getStatistic(xstat);
    Real stat(0);
    int nQuad = static_cast<int>(wts_.size());
    for (int i = 0; i < nQuad; ++i) {
      stat += wts_[i] * xstat[i];
    }
    return stat;
  }

  void reset(Teuchos::RCP<Vector<Real> > &x0, const Vector<Real> &x) {
    mqq_->reset(x0,x);
  }

  void reset(Teuchos::RCP<Vector<Real> > &x0, const Vector<Real> &x,
             Teuchos::RCP<Vector<Real> > &v0, const Vector<Real> &v) {
    mqq_->reset(x0,x,v0,v);
  }

  void update(const Real val, const Real weight) {
    mqq_->update(val,weight);
  }

  void update(const Real val, const Vector<Real> &g, const Real weight) {
    mqq_->update(val,g,weight);
  }

  void update(const Real val, const Vector<Real> &g, const Real gv, const Vector<Real> &hv,
              const Real weight) {
    mqq_->update(val,g,gv,hv,weight);
  }

  Real getValue(SampleGenerator<Real> &sampler) {
    return mqq_->getValue(sampler);
  }

  void getGradient(Vector<Real> &g, SampleGenerator<Real> &sampler) {
    mqq_->getGradient(g,sampler);
  }

  void getHessVec(Vector<Real> &hv, SampleGenerator<Real> &sampler) {
    mqq_->getHessVec(hv,sampler);
  }
};

}

#endif