44 #ifndef ROL_MIXEDQUANTILEQUADRANGLE_HPP
45 #define ROL_MIXEDQUANTILEQUADRANGLE_HPP
50 #include "ROL_ParameterList.hpp"
117 ROL_TEST_FOR_EXCEPTION((pSize!=cSize),std::invalid_argument,
118 ">>> ERROR (ROL::MixedCVaR): Probability and coefficient arrays have different sizes!");
119 Real sum(0),
zero(0), one(1);
120 for (
int i = 0; i < pSize; i++) {
121 ROL_TEST_FOR_EXCEPTION((
prob_[i]>one ||
prob_[i]<
zero), std::invalid_argument,
122 ">>> ERROR (ROL::MixedCVaR): Element of probability array out of range!");
123 ROL_TEST_FOR_EXCEPTION((
coeff_[i]>one ||
coeff_[i]<
zero), std::invalid_argument,
124 ">>> ERROR (ROL::MixedCVaR): Element of coefficient array out of range!");
127 ROL_TEST_FOR_EXCEPTION((std::abs(sum-one) > std::sqrt(ROL_EPSILON<Real>())),std::invalid_argument,
128 ">>> ERROR (ROL::MixedCVaR): Coefficients do not sum to one!");
129 ROL_TEST_FOR_EXCEPTION(
plusFunction_ == ROL::nullPtr, std::invalid_argument,
130 ">>> ERROR (ROL::MixedCVaR): PlusFunction pointer is null!");
138 ROL::ParameterList &list
139 = parlist.sublist(
"SOL").sublist(
"Risk Measure").sublist(
"Mixed CVaR");
141 prob_ = ROL::getArrayFromStringParameter<Real>(list,
"Probability Array");
142 coeff_ = ROL::getArrayFromStringParameter<Real>(list,
"Coefficient Array");
149 const std::vector<Real> &coeff,
162 if (xstat != nullPtr) {
163 for (
int i = 0; i <
size_; ++i) {
164 stat =
coeff_[i]*(*xstat)[i];
172 const std::vector<Real> &xstat,
176 for (
int i = 0; i <
size_; i++) {
183 const std::vector<Real> &xstat,
187 for (
int i = 0; i <
size_; i++) {
188 cvar +=
coeff_[i]*xstat[i];
195 const std::vector<Real> &xstat,
197 Real pf(0), c(0), one(1);
199 for (
int i = 0; i <
size_; i++) {
202 if (std::abs(c) >= ROL_EPSILON<Real>()) {
211 std::vector<Real> &gstat,
213 const std::vector<Real> &xstat,
216 for (
int i = 0; i <
size_; i++) {
224 const std::vector<Real> &vstat,
226 const std::vector<Real> &xstat,
228 Real pf1(0), pf2(0), c(0), one(1);
230 for (
int i = 0; i <
size_; i++) {
233 if (std::abs(pf2) >= ROL_EPSILON<Real>()) {
239 if (std::abs(pf1) >= ROL_EPSILON<Real>()) {
248 std::vector<Real> &hvstat,
250 const std::vector<Real> &vstat,
252 const std::vector<Real> &xstat,
Provides the interface to evaluate objective functions.
MixedCVaR(ROL::ParameterList &parlist)
void computeHessVec(Vector< Real > &hv, Objective< Real > &obj, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Real computeValue(Objective< Real > &obj, const Vector< Real > &x, Real &tol)
Ptr< Vector< Real > > hv_
void updateGradient(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal risk measure storage for gradient computation.
ROL::Ptr< PlusFunction< Real > > plusFunction_
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)
Update internal risk measure storage for Hessian-time-a-vector computation.
Ptr< Vector< Real > > dualVector_
std::vector< Real > coeff_
Defines the linear algebra or vector space interface.
void sumAll(Real *input, Real *output, int dim) const
Objective_SerialSimOpt(const Ptr< Obj > &obj, const V &ui) z0_ zero()
Real getValue(const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure value.
Real computeStatistic(const Ptr< std::vector< Real >> &xstat) const
Provides an interface for a convex combination of conditional value-at-risks.
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)
Return risk measure Hessian-times-a-vector.
void initialize(const Vector< Real > &x)
Initialize temporary variables.
MixedCVaR(const std::vector< Real > &prob, const std::vector< Real > &coeff, const ROL::Ptr< PlusFunction< Real > > &pf)
void computeGradient(Vector< Real > &g, Objective< Real > &obj, const Vector< Real > &x, Real &tol)
Real computeGradVec(Vector< Real > &g, Objective< Real > &obj, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Provides the interface to implement any functional that maps a random variable to a (extended) real n...
void initializeMCVAR(void)
std::vector< Real > prob_
void updateValue(Objective< Real > &obj, const Vector< Real > &x, const std::vector< Real > &xstat, Real &tol)
Update internal storage for value computation.
void getGradient(Vector< Real > &g, std::vector< Real > &gstat, const Vector< Real > &x, const std::vector< Real > &xstat, SampleGenerator< Real > &sampler)
Return risk measure (sub)gradient.
virtual void initialize(const Vector< Real > &x)
Initialize temporary variables.