10#ifndef ROL_QUANTILERADIUSQUADRANGLE_HPP
11#define ROL_QUANTILERADIUSQUADRANGLE_HPP
16#include "ROL_ParameterList.hpp"
51 ROL_TEST_FOR_EXCEPTION((
prob_>one ||
prob_<
zero), std::invalid_argument,
52 ">>> ERROR (ROL::QuantileRadius): Confidence level out of range!");
53 ROL_TEST_FOR_EXCEPTION((
coeff_<
zero), std::invalid_argument,
54 ">>> ERROR (ROL::QuantileRadius): Coefficient is negative!");
62 ROL::ParameterList &list
63 = parlist.sublist(
"SOL").sublist(
"Risk Measure").sublist(
"Quantile Radius");
65 prob_ = list.get<Real>(
"Confidence Level");
66 coeff_ = list.get<Real>(
"Coefficient");
80 vec_.assign(2,
static_cast<Real
>(0));
84 Real stat(0), half(0.5);
85 if (xstat != nullPtr) {
86 stat = half*((*xstat)[0] + (*xstat)[1]);
93 const std::vector<Real> &xstat,
95 const Real half(0.5), one(1);
96 Real val = computeValue(obj,x,tol);
103 const std::vector<Real> &xstat,
105 const Real half(0.5);
107 sampler.
sumAll(&val_,&cvar,1);
108 cvar += half*
coeff_*(xstat[0] + xstat[1]);
114 const std::vector<Real> &xstat,
115 Real &tol)
override {
116 const Real half(0.5), one(1);
117 Real val = computeValue(obj,x,tol);
123 computeGradient(*dualVector_,obj,x,tol);
124 g_->axpy(weight_ + c * (pf1 - pf2),*dualVector_);
128 std::vector<Real> &gstat,
130 const std::vector<Real> &xstat,
132 const Real half(0.5);
141 const std::vector<Real> &vstat,
143 const std::vector<Real> &xstat,
144 Real &tol)
override {
145 const Real half(0.5), one(1);
146 Real val = computeValue(obj,x,tol);
152 Real gv = computeGradVec(*dualVector_,obj,v,x,tol);
153 vec_[0] -= c*pf12*(gv-vstat[0]);
154 vec_[1] += c*pf22*(gv+vstat[1]);
155 hv_->axpy(c*(pf12*(gv-vstat[0]) + pf22*(gv+vstat[1])),*dualVector_);
156 computeHessVec(*dualVector_,obj,v,x,tol);
157 hv_->axpy(weight_ + c * (pf11 - pf21),*dualVector_);
161 std::vector<Real> &hvstat,
163 const std::vector<Real> &vstat,
165 const std::vector<Real> &xstat,
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0 zero)()
Provides the interface to evaluate objective functions.
void getGradient(Vector< Real > &g, std::vector< Real > &gstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler) override
void getHessVec(Vector< Real > &hv, std::vector< Real > &hvstat, const Vector< Real > &v, const std::vector< Real > &vstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler) override
void updateValue(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol) override
Real getValue(const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler) override
QuantileRadius(ROL::ParameterList &parlist)
void updateGradient(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol) override
QuantileRadius(const Real prob, const Real coeff, const Ptr< PlusFunction< Real > > &pf)
void updateHessVec(Objective< Real > &obj, const Vector< Real > &v, const std::vector< Real > &vstat, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol) override
Real computeStatistic(const Ptr< const std::vector< Real > > &xstat) const override
void initialize(const Vector< Real > &x) override
Ptr< PlusFunction< Real > > plusFunction_
Provides the interface to implement any functional that maps a random variable to a (extended) real n...
void sumAll(Real *input, Real *output, int dim) const
Defines the linear algebra or vector space interface.