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
48#include "ROL_PlusFunction.hpp"
49#include "ROL_RiskVector.hpp"
51
52namespace ROL {
53
54template <class Real>
55class SimulatedObjectiveCVaR : public Objective<Real> {
56private:
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
62public:
63
65
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 void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
77 pobj_->update(x,type,iter);
78 }
79
80 Real value(const Vector<Real> &x,
81 Real &tol) {
82 const Vector_SimOpt<Real> &uz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
83 ROL::Ptr<const Vector<Real> > uptr = uz.get_1();
84 ROL::Ptr<const Vector<Real> > zptr = uz.get_2();
85 const SimulatedVector<Real> &pu = dynamic_cast<const SimulatedVector<Real>&>(*uptr);
86 const RiskVector<Real> &rz = dynamic_cast<const RiskVector<Real>&>(*zptr);
87 Real t = (*rz.getStatistic(0))[0];
88 ROL::Ptr<const Vector<Real> > z = rz.getVector();
89
90 std::vector<Real> param;
91 Real weight(0), one(1);
92 Real val = 0;
93 Real tmpval = 0;
94 Real tmpsum = 0;
95 Real tmpplus = 0;
96 for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pu.numVectors(); ++i) {
97 param = sampler_->getMyPoint(static_cast<int>(i));
98 weight = sampler_->getMyWeight(static_cast<int>(i));
99 pobj_->setParameter(param);
100 //tmpval = pobj_->value(*(pu.get(i)), *zptr, tol);
101 pobj_->update(*(pu.get(i)), *z);
102 tmpval = pobj_->value(*(pu.get(i)), *z, tol);
103 tmpplus = pfunc_->evaluate(tmpval-t, 0);
104 tmpsum += tmpplus*weight;
105 }
106 sampler_->sumAll(&tmpsum, &val, 1);
107 val *= (one/(one-alpha_));
108 val += t;
109 return val;
110 }
111
112 virtual void gradient(Vector<Real> &g,
113 const Vector<Real> &x,
114 Real &tol) {
115 g.zero();
116 // split x
117 const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
118 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
119 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
120 const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
121 const RiskVector<Real> &rxz = dynamic_cast<const RiskVector<Real>&>(*xzptr);
122 Real xt = (*rxz.getStatistic(0))[0];
123 ROL::Ptr<const Vector<Real> > xz = rxz.getVector();
124 // split g
125 Vector_SimOpt<Real> &guz = dynamic_cast<Vector_SimOpt<Real>&>(g);
126 ROL::Ptr<Vector<Real> > guptr = guz.get_1();
127 ROL::Ptr<Vector<Real> > gzptr = guz.get_2();
128 SimulatedVector<Real> &pgu = dynamic_cast<SimulatedVector<Real>&>(*guptr);
129 RiskVector<Real> &rgz = dynamic_cast<RiskVector<Real>&>(*gzptr);
130 ROL::Ptr<Vector<Real> > gz = rgz.getVector();
131
132 std::vector<Real> param;
133 Real weight(0), one(1), sum(0), tmpsum(0), tmpval(0), tmpplus(0);
134 //ROL::Ptr<Vector<Real> > tmp1 = gzptr->clone();
135 //ROL::Ptr<Vector<Real> > tmp2 = gzptr->clone();
136 ROL::Ptr<Vector<Real> > tmp1 = gz->clone();
137 ROL::Ptr<Vector<Real> > tmp2 = gz->clone();
138 for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<pgu.numVectors(); ++i) {
139 param = sampler_->getMyPoint(static_cast<int>(i));
140 weight = sampler_->getMyWeight(static_cast<int>(i));
141 pobj_->setParameter(param);
142 pobj_->update(*(pxu.get(i)), *xz);
143 //tmpval = pobj_->value(*(pxu.get(i)), *xzptr, tol);
144 tmpval = pobj_->value(*(pxu.get(i)), *xz, tol);
145 tmpplus = pfunc_->evaluate(tmpval-xt, 1);
146 tmpsum += weight*tmpplus;
147 //Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
148 Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xz));
149 Vector_SimOpt<Real> gi(pgu.get(i), tmp1);
150 pobj_->gradient(gi, xi, tol);
151 gi.scale(weight*tmpplus);
152 tmp2->plus(*tmp1);
153 pgu.get(i)->scale(one/(one-alpha_));
154 }
155 //sampler_->sumAll(*tmp2, *gzptr);
156 //gzptr->scale(one/(one-alpha_));
157 sampler_->sumAll(*tmp2, *gz);
158 gz->scale(one/(one-alpha_));
159 sampler_->sumAll(&tmpsum, &sum, 1);
160 rgz.setStatistic(one - (one/(one-alpha_))*sum,0);
161 }
162
163/*
164 virtual void hessVec(Vector<Real> &hv,
165 const Vector<Real> &v,
166 const Vector<Real> &x,
167 Real &tol) {
168 hv.zero();
169 // split x
170 const Vector_SimOpt<Real> &xuz = dynamic_cast<const Vector_SimOpt<Real>&>(x);
171 ROL::Ptr<const Vector<Real> > xuptr = xuz.get_1();
172 ROL::Ptr<const Vector<Real> > xzptr = xuz.get_2();
173 const SimulatedVector<Real> &pxu = dynamic_cast<const SimulatedVector<Real>&>(*xuptr);
174 // split v
175 const Vector_SimOpt<Real> &vuz = dynamic_cast<const Vector_SimOpt<Real>&>(v);
176 ROL::Ptr<const Vector<Real> > vuptr = vuz.get_1();
177 ROL::Ptr<const Vector<Real> > vzptr = vuz.get_2();
178 const SimulatedVector<Real> &pvu = dynamic_cast<const SimulatedVector<Real>&>(*vuptr);
179 // split hv
180 Vector_SimOpt<Real> &hvuz = dynamic_cast<Vector_SimOpt<Real>&>(hv);
181 ROL::Ptr<Vector<Real> > hvuptr = hvuz.get_1();
182 ROL::Ptr<Vector<Real> > hvzptr = hvuz.get_2();
183 SimulatedVector<Real> &phvu = dynamic_cast<SimulatedVector<Real>&>(*hvuptr);
184
185 std::vector<Real> param;
186 Real weight(0);
187 ROL::Ptr<Vector<Real> > tmp1 = hvzptr->clone();
188 ROL::Ptr<Vector<Real> > tmp2 = hvzptr->clone();
189 for (typename std::vector<SimulatedVector<Real> >::size_type i=0; i<phvu.numVectors(); ++i) {
190 param = sampler_->getMyPoint(static_cast<int>(i));
191 weight = sampler_->getMyWeight(static_cast<int>(i));
192 pobj_->setParameter(param);
193 Vector_SimOpt<Real> xi(ROL::constPtrCast<Vector<Real> >(pxu.get(i)), ROL::constPtrCast<Vector<Real> >(xzptr));
194 Vector_SimOpt<Real> vi(ROL::constPtrCast<Vector<Real> >(pvu.get(i)), ROL::constPtrCast<Vector<Real> >(vzptr));
195 Vector_SimOpt<Real> hvi(phvu.get(i), tmp1);
196 pobj_->update(xi);
197 pobj_->hessVec(hvi, vi, xi, tol);
198 hvi.scale(weight);
199 tmp2->plus(*tmp1);
200 }
201 sampler_->sumAll(*tmp2, *hvzptr);
202 }
203*/
204
205}; // class SimulatedObjective
206
207} // namespace ROL
208
209#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.