ROL
ROL_SimulatedObjectiveCVaR.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Rapid Optimization Library (ROL) Package
5 // Copyright (2014) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact lead developers:
38 // Drew Kouri (dpkouri@sandia.gov) and
39 // Denis Ridzal (dridzal@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
44 #ifndef ROL_SIMULATED_OBJECTIVE_CVAR_H
45 #define ROL_SIMULATED_OBJECTIVE_CVAR_H
46 
47 #include "ROL_SimulatedVector.hpp"
48 #include "ROL_PlusFunction.hpp"
49 #include "ROL_RiskVector.hpp"
50 #include "ROL_Objective_SimOpt.hpp"
51 
52 namespace ROL {
53 
54 template <class Real>
55 class SimulatedObjectiveCVaR : public Objective<Real> {
56 private:
57  const ROL::Ptr<SampleGenerator<Real> > sampler_;
58  const ROL::Ptr<Objective_SimOpt<Real> > pobj_;
59  const ROL::Ptr<PlusFunction<Real> > pfunc_;
60  const Real alpha_;
61 
62 public:
63 
65 
66  SimulatedObjectiveCVaR(const ROL::Ptr<SampleGenerator<Real> > & sampler,
67  const ROL::Ptr<Objective_SimOpt<Real> > & pobj,
68  const ROL::Ptr<PlusFunction<Real> > & pfunc,
69  const Real & alpha)
70  : sampler_(sampler), pobj_(pobj), pfunc_(pfunc), alpha_(alpha) {}
71 
72  void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
73  pobj_->update(x,flag,iter);
74  }
75 
76  Real value(const Vector<Real> &x,
77  Real &tol) {
78  const Vector_SimOpt<Real> &uz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
79  ROL::Ptr<const Vector<Real> > uptr = uz.get_1();
80  ROL::Ptr<const Vector<Real> > zptr = uz.get_2();
81  const SimulatedVector<Real> &pu = dynamic_cast<const SimulatedVector<Real>&>(*uptr);
82  const RiskVector<Real> &rz = dynamic_cast<const RiskVector<Real>&>(*zptr);
83  Real t = (*rz.getStatistic(0))[0];
84  ROL::Ptr<const Vector<Real> > z = rz.getVector();
85 
86  std::vector<Real> param;
87  Real weight(0), one(1);
88  Real val = 0;
89  Real tmpval = 0;
90  Real tmpsum = 0;
91  Real tmpplus = 0;
92  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pu.numVectors(); ++i) {
93  param = sampler_->getMyPoint(static_cast<int>(i));
94  weight = sampler_->getMyWeight(static_cast<int>(i));
95  pobj_->setParameter(param);
96  //tmpval = pobj_->value(*(pu.get(i)), *zptr, tol);
97  pobj_->update(*(pu.get(i)), *z);
98  tmpval = pobj_->value(*(pu.get(i)), *z, tol);
99  tmpplus = pfunc_->evaluate(tmpval-t, 0);
100  tmpsum += tmpplus*weight;
101  }
102  sampler_->sumAll(&tmpsum, &val, 1);
103  val *= (one/(one-alpha_));
104  val += t;
105  return val;
106  }
107 
108  virtual void gradient(Vector<Real> &g,
109  const Vector<Real> &x,
110  Real &tol) {
111  g.zero();
112  // split x
113  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
114  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
115  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
116  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
117  const RiskVector<Real> &rxz = dynamic_cast<const RiskVector<Real>&>(*xzptr);
118  Real xt = (*rxz.getStatistic(0))[0];
119  ROL::Ptr<const Vector<Real> > xz = rxz.getVector();
120  // split g
121  Vector_SimOpt<Real> &guz = dynamic_cast<Vector_SimOpt<Real>&>(g);
122  ROL::Ptr<Vector<Real> > guptr = guz.get_1();
123  ROL::Ptr<Vector<Real> > gzptr = guz.get_2();
124  SimulatedVector<Real> &pgu = dynamic_cast<SimulatedVector<Real>&>(*guptr);
125  RiskVector<Real> &rgz = dynamic_cast<RiskVector<Real>&>(*gzptr);
126  ROL::Ptr<Vector<Real> > gz = rgz.getVector();
127 
128  std::vector<Real> param;
129  Real weight(0), one(1), sum(0), tmpsum(0), tmpval(0), tmpplus(0);
130  //ROL::Ptr<Vector<Real> > tmp1 = gzptr->clone();
131  //ROL::Ptr<Vector<Real> > tmp2 = gzptr->clone();
132  ROL::Ptr<Vector<Real> > tmp1 = gz->clone();
133  ROL::Ptr<Vector<Real> > tmp2 = gz->clone();
134  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pgu.numVectors(); ++i) {
135  param = sampler_->getMyPoint(static_cast<int>(i));
136  weight = sampler_->getMyWeight(static_cast<int>(i));
137  pobj_->setParameter(param);
138  pobj_->update(*(pxu.get(i)), *xz);
139  //tmpval = pobj_->value(*(pxu.get(i)), *xzptr, tol);
140  tmpval = pobj_->value(*(pxu.get(i)), *xz, tol);
141  tmpplus = pfunc_->evaluate(tmpval-xt, 1);
142  tmpsum += weight*tmpplus;
143  //Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
144  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xz));
145  Vector_SimOpt<Real> gi(pgu.get(i), tmp1);
146  pobj_->gradient(gi, xi, tol);
147  gi.scale(weight*tmpplus);
148  tmp2->plus(*tmp1);
149  pgu.get(i)->scale(one/(one-alpha_));
150  }
151  //sampler_->sumAll(*tmp2, *gzptr);
152  //gzptr->scale(one/(one-alpha_));
153  sampler_->sumAll(*tmp2, *gz);
154  gz->scale(one/(one-alpha_));
155  sampler_->sumAll(&tmpsum, &sum, 1);
156  rgz.setStatistic(one - (one/(one-alpha_))*sum,0);
157  }
158 
159 /*
160  virtual void hessVec(Vector<Real> &hv,
161  const Vector<Real> &v,
162  const Vector<Real> &x,
163  Real &tol) {
164  hv.zero();
165  // split x
166  const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
167  ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
168  ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
169  const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
170  // split v
171  const Vector_SimOpt<Real> &vuz = dynamic_cast<const Vector_SimOpt<Real>&>(v);
172  ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
173  ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
174  const SimulatedVector<Real> &pvu = dynamic_cast<const SimulatedVector<Real>&>(*vuptr);
175  // split hv
176  Vector_SimOpt<Real> &hvuz = dynamic_cast<Vector_SimOpt<Real>&>(hv);
177  ROL::Ptr<Vector<Real> > hvuptr = hvuz.get_1();
178  ROL::Ptr<Vector<Real> > hvzptr = hvuz.get_2();
179  SimulatedVector<Real> &phvu = dynamic_cast<SimulatedVector<Real>&>(*hvuptr);
180 
181  std::vector<Real> param;
182  Real weight(0);
183  ROL::Ptr<Vector<Real> > tmp1 = hvzptr->clone();
184  ROL::Ptr<Vector<Real> > tmp2 = hvzptr->clone();
185  for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<phvu.numVectors(); ++i) {
186  param = sampler_->getMyPoint(static_cast<int>(i));
187  weight = sampler_->getMyWeight(static_cast<int>(i));
188  pobj_->setParameter(param);
189  Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
190  Vector_SimOpt<Real> vi(ROL::constPtrCast<Vector<Real> >(pvu.get(i)), ROL::constPtrCast<Vector<Real> >(vzptr));
191  Vector_SimOpt<Real> hvi(phvu.get(i), tmp1);
192  pobj_->update(xi);
193  pobj_->hessVec(hvi, vi, xi, tol);
194  hvi.scale(weight);
195  tmp2->plus(*tmp1);
196  }
197  sampler_->sumAll(*tmp2, *hvzptr);
198  }
199 */
200 
201 }; // class SimulatedObjective
202 
203 } // namespace ROL
204 
205 #endif
Provides the interface to evaluate objective functions.
Provides the interface to evaluate simulation-based objective functions.
ROL::Ptr< std::vector< Real > > getStatistic(const int comp=0, const int index=0)
typename PV< Real >::size_type size_type
const ROL::Ptr< PlusFunction< Real > > pfunc_
ROL::Ptr< const Vector< Real > > get_2() const
Defines the linear algebra or vector space interface for simulation-based optimization.
ROL::Ptr< const Vector< Real > > getVector(void) const
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
SimulatedObjectiveCVaR(const ROL::Ptr< SampleGenerator< Real > > &sampler, const ROL::Ptr< Objective_SimOpt< Real > > &pobj, const ROL::Ptr< PlusFunction< Real > > &pfunc, const Real &alpha)
void setStatistic(const Real stat, const int comp=0, const int index=0)
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:167
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:80
Defines the linear algebra of a vector space on a generic partitioned vector where the individual vec...
ROL::Ptr< const Vector< Real > > get(size_type i) const
size_type numVectors() const
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
const ROL::Ptr< SampleGenerator< Real > > sampler_
Real value(const Vector< Real > &x, Real &tol)
Compute value.
const ROL::Ptr< Objective_SimOpt< Real > > pobj_
ROL::Ptr< const Vector< Real > > get_1() const