ROL
ROL_SimulatedObjectiveCVaR.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_SIMULATED_OBJECTIVE_CVAR_H
11#define ROL_SIMULATED_OBJECTIVE_CVAR_H
12
14#include "ROL_PlusFunction.hpp"
15#include "ROL_RiskVector.hpp"
17
18namespace ROL {
19
20template <class Real>
21class SimulatedObjectiveCVaR : public Objective<Real> {
22private:
23 const ROL::Ptr<SampleGenerator<Real> > sampler_;
24 const ROL::Ptr<Objective_SimOpt<Real> > pobj_;
25 const ROL::Ptr<PlusFunction<Real> > pfunc_;
26 const Real alpha_;
27
28public:
29
31
33 const ROL::Ptr<Objective_SimOpt<Real> > & pobj,
34 const ROL::Ptr<PlusFunction<Real> > & pfunc,
35 const Real & alpha)
36 : sampler_(sampler), pobj_(pobj), pfunc_(pfunc), alpha_(alpha) {}
37
38 void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
39 pobj_->update(x,flag,iter);
40 }
41
42 void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
43 pobj_->update(x,type,iter);
44 }
45
46 Real value(const Vector<Real> &x,
47 Real &tol) {
48 const Vector_SimOpt<Real> &uz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
49 ROL::Ptr<const Vector<Real> > uptr = uz.get_1();
50 ROL::Ptr<const Vector<Real> > zptr = uz.get_2();
51 const SimulatedVector<Real> &pu = dynamic_cast<const SimulatedVector<Real>&>(*uptr);
52 const RiskVector<Real> &rz = dynamic_cast<const RiskVector<Real>&>(*zptr);
53 Real t = (*rz.getStatistic(0))[0];
54 ROL::Ptr<const Vector<Real> > z = rz.getVector();
55
56 std::vector<Real> param;
57 Real weight(0), one(1);
58 Real val = 0;
59 Real tmpval = 0;
60 Real tmpsum = 0;
61 Real tmpplus = 0;
62 for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pu.numVectors(); ++i) {
63 param = sampler_->getMyPoint(static_cast<int>(i));
64 weight = sampler_->getMyWeight(static_cast<int>(i));
65 pobj_->setParameter(param);
66 //tmpval = pobj_->value(*(pu.get(i)), *zptr, tol);
67 pobj_->update(*(pu.get(i)), *z);
68 tmpval = pobj_->value(*(pu.get(i)), *z, tol);
69 tmpplus = pfunc_->evaluate(tmpval-t, 0);
70 tmpsum += tmpplus*weight;
71 }
72 sampler_->sumAll(&tmpsum, &val, 1);
73 val *= (one/(one-alpha_));
74 val += t;
75 return val;
76 }
77
78 virtual void gradient(Vector<Real> &g,
79 const Vector<Real> &x,
80 Real &tol) {
81 g.zero();
82 // split x
83 const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
84 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
85 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
86 const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
87 const RiskVector<Real> &rxz = dynamic_cast<const RiskVector<Real>&>(*xzptr);
88 Real xt = (*rxz.getStatistic(0))[0];
89 ROL::Ptr<const Vector<Real> > xz = rxz.getVector();
90 // split g
91 Vector_SimOpt<Real> &guz = dynamic_cast<Vector_SimOpt<Real>&>(g);
92 ROL::Ptr<Vector<Real> > guptr = guz.get_1();
93 ROL::Ptr<Vector<Real> > gzptr = guz.get_2();
94 SimulatedVector<Real> &pgu = dynamic_cast<SimulatedVector<Real>&>(*guptr);
95 RiskVector<Real> &rgz = dynamic_cast<RiskVector<Real>&>(*gzptr);
96 ROL::Ptr<Vector<Real> > gz = rgz.getVector();
97
98 std::vector<Real> param;
99 Real weight(0), one(1), sum(0), tmpsum(0), tmpval(0), tmpplus(0);
100 //ROL::Ptr<Vector<Real> > tmp1 = gzptr->clone();
101 //ROL::Ptr<Vector<Real> > tmp2 = gzptr->clone();
102 ROL::Ptr<Vector<Real> > tmp1 = gz->clone();
103 ROL::Ptr<Vector<Real> > tmp2 = gz->clone();
104 for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pgu.numVectors(); ++i) {
105 param = sampler_->getMyPoint(static_cast<int>(i));
106 weight = sampler_->getMyWeight(static_cast<int>(i));
107 pobj_->setParameter(param);
108 pobj_->update(*(pxu.get(i)), *xz);
109 //tmpval = pobj_->value(*(pxu.get(i)), *xzptr, tol);
110 tmpval = pobj_->value(*(pxu.get(i)), *xz, tol);
111 tmpplus = pfunc_->evaluate(tmpval-xt, 1);
112 tmpsum += weight*tmpplus;
113 //Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
114 Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xz));
115 Vector_SimOpt<Real> gi(pgu.get(i), tmp1);
116 pobj_->gradient(gi, xi, tol);
117 gi.scale(weight*tmpplus);
118 tmp2->plus(*tmp1);
119 pgu.get(i)->scale(one/(one-alpha_));
120 }
121 //sampler_->sumAll(*tmp2, *gzptr);
122 //gzptr->scale(one/(one-alpha_));
123 sampler_->sumAll(*tmp2, *gz);
124 gz->scale(one/(one-alpha_));
125 sampler_->sumAll(&tmpsum, &sum, 1);
126 rgz.setStatistic(one - (one/(one-alpha_))*sum,0);
127 }
128
129/*
130 virtual void hessVec(Vector<Real> &hv,
131 const Vector<Real> &v,
132 const Vector<Real> &x,
133 Real &tol) {
134 hv.zero();
135 // split x
136 const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
137 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
138 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
139 const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
140 // split v
141 const Vector_SimOpt<Real> &vuz = dynamic_cast<const Vector_SimOpt<Real>&>(v);
142 ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
143 ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
144 const SimulatedVector<Real> &pvu = dynamic_cast<const SimulatedVector<Real>&>(*vuptr);
145 // split hv
146 Vector_SimOpt<Real> &hvuz = dynamic_cast<Vector_SimOpt<Real>&>(hv);
147 ROL::Ptr<Vector<Real> > hvuptr = hvuz.get_1();
148 ROL::Ptr<Vector<Real> > hvzptr = hvuz.get_2();
149 SimulatedVector<Real> &phvu = dynamic_cast<SimulatedVector<Real>&>(*hvuptr);
150
151 std::vector<Real> param;
152 Real weight(0);
153 ROL::Ptr<Vector<Real> > tmp1 = hvzptr->clone();
154 ROL::Ptr<Vector<Real> > tmp2 = hvzptr->clone();
155 for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<phvu.numVectors(); ++i) {
156 param = sampler_->getMyPoint(static_cast<int>(i));
157 weight = sampler_->getMyWeight(static_cast<int>(i));
158 pobj_->setParameter(param);
159 Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
160 Vector_SimOpt<Real> vi(ROL::constPtrCast<Vector<Real> >(pvu.get(i)), ROL::constPtrCast<Vector<Real> >(vzptr));
161 Vector_SimOpt<Real> hvi(phvu.get(i), tmp1);
162 pobj_->update(xi);
163 pobj_->hessVec(hvi, vi, xi, tol);
164 hvi.scale(weight);
165 tmp2->plus(*tmp1);
166 }
167 sampler_->sumAll(*tmp2, *hvzptr);
168 }
169*/
170
171}; // class SimulatedObjective
172
173} // namespace ROL
174
175#endif
typename PV< Real >::size_type size_type
Provides the interface to evaluate simulation-based objective functions.
Provides the interface to evaluate objective functions.
Ptr< std::vector< Real > > getStatistic(const int comp=0, const int index=0)
void setStatistic(const Real stat, const int comp=0, const int index=0)
Ptr< const Vector< Real > > getVector(void) const
Real value(const Vector< Real > &x, Real &tol)
Compute value.
const ROL::Ptr< PlusFunction< Real > > pfunc_
virtual void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
SimulatedObjectiveCVaR(const ROL::Ptr< SampleGenerator< Real > > &sampler, const ROL::Ptr< Objective_SimOpt< Real > > &pobj, const ROL::Ptr< PlusFunction< Real > > &pfunc, const Real &alpha)
const ROL::Ptr< Objective_SimOpt< Real > > pobj_
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
const ROL::Ptr< SampleGenerator< Real > > sampler_
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
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
Defines the linear algebra or vector space interface for simulation-based optimization.
ROL::Ptr< const Vector< Real > > get_2() const
ROL::Ptr< const Vector< Real > > get_1() const
void scale(const Real alpha)
Compute where .
Defines the linear algebra or vector space interface.
virtual void zero()
Set to zero vector.