ROL
ROL_MixedCVaR.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_MIXEDQUANTILEQUADRANGLE_HPP
11 #define ROL_MIXEDQUANTILEQUADRANGLE_HPP
12 
14 #include "ROL_PlusFunction.hpp"
15 
16 #include "ROL_ParameterList.hpp"
17 
51 namespace ROL {
52 
53 template<class Real>
54 class MixedCVaR : public RandVarFunctional<Real> {
55 private:
56  ROL::Ptr<PlusFunction<Real> > plusFunction_;
57  std::vector<Real> prob_;
58  std::vector<Real> coeff_;
59  std::vector<Real> vec_;
60  int size_;
61 
67 
70 
75 
76  void initializeMCVAR(void) {
77  size_ = prob_.size();
78  vec_.clear(); vec_.resize(size_,static_cast<Real>(0));
79  }
80 
81  void checkInputs(void) {
82  int pSize = prob_.size(), cSize = coeff_.size();
83  ROL_TEST_FOR_EXCEPTION((pSize!=cSize),std::invalid_argument,
84  ">>> ERROR (ROL::MixedCVaR): Probability and coefficient arrays have different sizes!");
85  Real sum(0), zero(0), one(1);
86  for (int i = 0; i < pSize; i++) {
87  ROL_TEST_FOR_EXCEPTION((prob_[i]>one || prob_[i]<zero), std::invalid_argument,
88  ">>> ERROR (ROL::MixedCVaR): Element of probability array out of range!");
89  ROL_TEST_FOR_EXCEPTION((coeff_[i]>one || coeff_[i]<zero), std::invalid_argument,
90  ">>> ERROR (ROL::MixedCVaR): Element of coefficient array out of range!");
91  sum += coeff_[i];
92  }
93  ROL_TEST_FOR_EXCEPTION((std::abs(sum-one) > std::sqrt(ROL_EPSILON<Real>())),std::invalid_argument,
94  ">>> ERROR (ROL::MixedCVaR): Coefficients do not sum to one!");
95  ROL_TEST_FOR_EXCEPTION(plusFunction_ == ROL::nullPtr, std::invalid_argument,
96  ">>> ERROR (ROL::MixedCVaR): PlusFunction pointer is null!");
98  }
99 
100 public:
101 
102  MixedCVaR( ROL::ParameterList &parlist )
103  : RandVarFunctional<Real>() {
104  ROL::ParameterList &list
105  = parlist.sublist("SOL").sublist("Risk Measure").sublist("Mixed CVaR");
106  // Grab probability and coefficient arrays
107  prob_ = ROL::getArrayFromStringParameter<Real>(list,"Probability Array");
108  coeff_ = ROL::getArrayFromStringParameter<Real>(list,"Coefficient Array");
109  plusFunction_ = ROL::makePtr<PlusFunction<Real>>(list);
110  // Check inputs
111  checkInputs();
112  }
113 
114  MixedCVaR(const std::vector<Real> &prob,
115  const std::vector<Real> &coeff,
116  const ROL::Ptr<PlusFunction<Real> > &pf )
117  : RandVarFunctional<Real>(), plusFunction_(pf), prob_(prob), coeff_(coeff) {
118  checkInputs();
119  }
120 
121  void initialize(const Vector<Real> &x) {
123  vec_.assign(size_,static_cast<Real>(0));
124  }
125 
126  Real computeStatistic(const Ptr<const std::vector<Real>> &xstat) const override {
127  Real stat(0);
128  if (xstat != nullPtr) {
129  for (int i = 0; i < size_; ++i) {
130  stat = coeff_[i]*(*xstat)[i];
131  }
132  }
133  return stat;
134  }
135 
137  const Vector<Real> &x,
138  const std::vector<Real> &xstat,
139  Real &tol) {
140  Real pf(0), one(1);
141  Real val = computeValue(obj,x,tol);
142  for (int i = 0; i < size_; i++) {
143  pf = plusFunction_->evaluate(val-xstat[i],0);
144  val_ += weight_*coeff_[i]/(one-prob_[i])*pf;
145  }
146  }
147 
148  Real getValue(const Vector<Real> &x,
149  const std::vector<Real> &xstat,
150  SampleGenerator<Real> &sampler) {
151  Real cvar(0);
152  sampler.sumAll(&val_,&cvar,1);
153  for (int i = 0; i < size_; i++) {
154  cvar += coeff_[i]*xstat[i];
155  }
156  return cvar;
157  }
158 
160  const Vector<Real> &x,
161  const std::vector<Real> &xstat,
162  Real &tol) {
163  Real pf(0), c(0), one(1);
164  Real val = computeValue(obj,x,tol);
165  for (int i = 0; i < size_; i++) {
166  pf = plusFunction_->evaluate(val-xstat[i],1);
167  c = weight_*coeff_[i]/(one-prob_[i])*pf;
168  if (std::abs(c) >= ROL_EPSILON<Real>()) {
169  vec_[i] -= c;
170  computeGradient(*dualVector_,obj,x,tol);
171  g_->axpy(c,*dualVector_);
172  }
173  }
174  }
175 
177  std::vector<Real> &gstat,
178  const Vector<Real> &x,
179  const std::vector<Real> &xstat,
180  SampleGenerator<Real> &sampler) {
181  sampler.sumAll(&vec_[0],&gstat[0],size_);
182  for (int i = 0; i < size_; i++) {
183  gstat[i] += coeff_[i];
184  }
185  sampler.sumAll(*g_,g);
186  }
187 
189  const Vector<Real> &v,
190  const std::vector<Real> &vstat,
191  const Vector<Real> &x,
192  const std::vector<Real> &xstat,
193  Real &tol) {
194  Real pf1(0), pf2(0), c(0), one(1);
195  Real val = computeValue(obj,x,tol);
196  for (int i = 0; i < size_; i++) {
197  pf1 = plusFunction_->evaluate(val-xstat[i],1);
198  pf2 = plusFunction_->evaluate(val-xstat[i],2);
199  if (std::abs(pf2) >= ROL_EPSILON<Real>()) {
200  Real gv = computeGradVec(*dualVector_,obj,v,x,tol);
201  c = weight_*coeff_[i]/(one-prob_[i])*pf2*(gv-vstat[i]);
202  vec_[i] -= c;
203  hv_->axpy(c,*dualVector_);
204  }
205  if (std::abs(pf1) >= ROL_EPSILON<Real>()) {
206  c = weight_*coeff_[i]/(one-prob_[i])*pf1;
207  computeHessVec(*dualVector_,obj,v,x,tol);
209  }
210  }
211  }
212 
214  std::vector<Real> &hvstat,
215  const Vector<Real> &v,
216  const std::vector<Real> &vstat,
217  const Vector<Real> &x,
218  const std::vector<Real> &xstat,
219  SampleGenerator<Real> &sampler) {
220  sampler.sumAll(&vec_[0],&hvstat[0],size_);
221  sampler.sumAll(*hv_,hv);
222  }
223 };
224 
225 }
226 
227 #endif
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)
Ptr< Vector< Real > > g_
Real computeValue(Objective< Real > &obj, const Vector< Real > &x, Real &tol)
std::vector< Real > vec_
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.
Definition: ROL_Vector.hpp:46
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.
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.
void checkInputs(void)
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.
Real computeStatistic(const Ptr< const std::vector< Real >> &xstat) const override
Compute statistic.