ROL
ROL_SecondOrderCVaR.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_SUPERQUANTILEQUADRANGLE_HPP
11#define ROL_SUPERQUANTILEQUADRANGLE_HPP
12
13#include "ROL_SpectralRisk.hpp"
14#include "ROL_GaussLegendreQuadrature.hpp"
15#include "ROL_Fejer2Quadrature.hpp"
16
48namespace ROL {
49
50template<class Real>
51class SecondOrderCVaR : public SpectralRisk<Real> {
52private:
53 ROL::Ptr<PlusFunction<Real> > plusFunction_;
54
55 Real alpha_;
56 int nQuad_;
58
59 std::vector<Real> wts_;
60 std::vector<Real> pts_;
61
62 void checkInputs(void) const {
63 ROL_TEST_FOR_EXCEPTION((alpha_ < 0 || alpha_ >= 1), std::invalid_argument,
64 ">>> ERROR (ROL::SecondOrderCVaR): Confidence level not between 0 and 1!");
65 ROL_TEST_FOR_EXCEPTION(plusFunction_ == ROL::nullPtr, std::invalid_argument,
66 ">>> ERROR (ROL::SecondOrderCVaR): PlusFunction pointer is null!");
67 }
68
69 void initializeQuad(void) {
70 ROL::Ptr<Quadrature1D<Real> > quad;
71 if ( useGauss_ ) {
72 quad = ROL::makePtr<GaussLegendreQuadrature<Real>>(nQuad_);
73 }
74 else {
75 quad = ROL::makePtr<Fejer2Quadrature<Real>>(nQuad_);
76 }
77 // quad->test();
78 quad->get(pts_,wts_);
79 Real sum(0), half(0.5), one(1);
80 for (int i = 0; i < nQuad_; ++i) {
81 sum += wts_[i];
82 }
83 for (int i = 0; i < nQuad_; ++i) {
84 wts_[i] /= sum;
85 pts_[i] = one - alpha_*(half*(pts_[i] + one));
86 }
88 }
89
90public:
91
92 SecondOrderCVaR( ROL::ParameterList &parlist )
93 : SpectralRisk<Real>() {
94 ROL::ParameterList &list
95 = parlist.sublist("SOL").sublist("Risk Measure").sublist("Second Order CVaR");
96 // Grab confidence level and quadrature order
97 alpha_ = list.get<Real>("Confidence Level");
98 nQuad_ = list.get("Number of Quadrature Points",5);
99 useGauss_ = list.get("Use Gauss-Legendre Quadrature",true);
100 plusFunction_ = ROL::makePtr<PlusFunction<Real>>(list);
101 // Check inputs
102 checkInputs();
104 }
105
106 SecondOrderCVaR(const Real alpha,
107 const int nQuad,
108 const ROL::Ptr<PlusFunction<Real> > &pf,
109 const bool useGauss = true)
110 : SpectralRisk<Real>(), plusFunction_(pf),
111 alpha_(alpha), nQuad_(nQuad), useGauss_(useGauss) {
112 // Check inputs
113 checkInputs();
115 }
116};
117
118}
119
120#endif
Provides an interface for the risk measure associated with the super quantile quadrangle.
SecondOrderCVaR(ROL::ParameterList &parlist)
SecondOrderCVaR(const Real alpha, const int nQuad, const ROL::Ptr< PlusFunction< Real > > &pf, const bool useGauss=true)
std::vector< Real > pts_
std::vector< Real > wts_
ROL::Ptr< PlusFunction< Real > > plusFunction_
Provides an interface for spectral risk measures.
void buildMixedQuantile(const std::vector< Real > &pts, const std::vector< Real > &wts, const ROL::Ptr< PlusFunction< Real > > &pf)