ROL
ROL_RiskNeutralObjective.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_RISKNEUTRALOBJECTIVE_HPP
45#define ROL_RISKNEUTRALOBJECTIVE_HPP
46
47#include "ROL_Vector.hpp"
48#include "ROL_Objective.hpp"
52
53namespace ROL {
54
55template<class Real>
56class RiskNeutralObjective : public Objective<Real> {
57private:
58 Ptr<Objective<Real>> ParametrizedObjective_;
59 Ptr<SampleGenerator<Real>> ValueSampler_;
60 Ptr<SampleGenerator<Real>> GradientSampler_;
61 Ptr<SampleGenerator<Real>> HessianSampler_;
62
63 Real value_;
64 Ptr<Vector<Real>> gradient_;
65 Ptr<Vector<Real>> pointDual_;
66 Ptr<Vector<Real>> sumDual_;
67
70
71 //std::map<std::vector<Real>,Real> value_storage_;
72 //std::map<std::vector<Real>,Ptr<Vector<Real>>> gradient_storage_;
73 Ptr<ScalarController<Real>> value_storage_;
74 Ptr<VectorController<Real>> gradient_storage_;
75
76 void initialize(const Vector<Real> &x) {
77 if ( firstUpdate_ ) {
78 gradient_ = (x.dual()).clone();
79 pointDual_ = (x.dual()).clone();
80 sumDual_ = (x.dual()).clone();
81 firstUpdate_ = false;
82 }
83 }
84
85 void getValue(Real &val, const Vector<Real> &x,
86 const std::vector<Real> &param, Real &tol) {
87 bool isComputed = false;
88 if ( storage_) {
89 isComputed = value_storage_->get(val,param);
90 }
91 if (!isComputed || !storage_) {
92 ParametrizedObjective_->setParameter(param);
93 val = ParametrizedObjective_->value(x,tol);
94 if ( storage_ ) {
95 value_storage_->set(val,param);
96 }
97 }
98 }
99
101 const std::vector<Real> &param, Real &tol) {
102 bool isComputed = false;
103 if ( storage_) {
104 isComputed = gradient_storage_->get(g,param);
105 }
106 if (!isComputed || !storage_) {
107 ParametrizedObjective_->setParameter(param);
108 ParametrizedObjective_->gradient(g,x,tol);
109 if ( storage_ ) {
110 gradient_storage_->set(g,param);
111 }
112 }
113 }
114
115 void getHessVec(Vector<Real> &hv, const Vector<Real> &v, const Vector<Real> &x,
116 const std::vector<Real> &param, Real &tol) {
117 ParametrizedObjective_->setParameter(param);
118 ParametrizedObjective_->hessVec(hv,v,x,tol);
119 }
120
121
122public:
124 const Ptr<SampleGenerator<Real>> &vsampler,
125 const Ptr<SampleGenerator<Real>> &gsampler,
126 const Ptr<SampleGenerator<Real>> &hsampler,
127 const bool storage = true )
129 ValueSampler_(vsampler), GradientSampler_(gsampler), HessianSampler_(hsampler),
130 firstUpdate_(true), storage_(storage) {
131 value_storage_ = makePtr<ScalarController<Real>>();
132 gradient_storage_ = makePtr<VectorController<Real>>();
133 }
134
136 const Ptr<SampleGenerator<Real>> &vsampler,
137 const Ptr<SampleGenerator<Real>> &gsampler,
138 const bool storage = true )
140 ValueSampler_(vsampler), GradientSampler_(gsampler), HessianSampler_(gsampler),
141 firstUpdate_(true), storage_(storage) {
142 value_storage_ = makePtr<ScalarController<Real>>();
143 gradient_storage_ = makePtr<VectorController<Real>>();
144 }
145
147 const Ptr<SampleGenerator<Real>> &sampler,
148 const bool storage = true )
150 ValueSampler_(sampler), GradientSampler_(sampler), HessianSampler_(sampler),
151 firstUpdate_(true), storage_(storage) {
152 value_storage_ = makePtr<ScalarController<Real>>();
153 gradient_storage_ = makePtr<VectorController<Real>>();
154 }
155
156 void update( const Vector<Real> &x, UpdateType type, int iter = -1 ) {
157 initialize(x);
158// ParametrizedObjective_->update(x,(flag && iter>=0),iter);
159 ParametrizedObjective_->update(x,type,iter);
160 ValueSampler_->update(x);
161 value_ = static_cast<Real>(0);
162 if ( storage_ ) {
163 value_storage_->objectiveUpdate(type);
164 gradient_storage_->objectiveUpdate(type);
165 }
166 if ( type != UpdateType::Trial && type != UpdateType::Revert ) { //&& iter>=0 ) {
167 GradientSampler_->update(x);
168 HessianSampler_->update(x);
169 gradient_->zero();
170 }
171 }
172
173 void update( const Vector<Real> &x, bool flag = true, int iter = -1 ) {
174 initialize(x);
175// ParametrizedObjective_->update(x,(flag && iter>=0),iter);
176 ParametrizedObjective_->update(x,flag,iter);
177 ValueSampler_->update(x);
178 value_ = static_cast<Real>(0);
179 if ( storage_ ) {
180 value_storage_->objectiveUpdate(true);
181 }
182 //if ( flag ) { //&& iter>=0 ) {
183 GradientSampler_->update(x);
184 HessianSampler_->update(x);
185 gradient_->zero();
186 if ( storage_ ) {
187 gradient_storage_->objectiveUpdate(true);
188 }
189 //}
190 }
191
192 Real value( const Vector<Real> &x, Real &tol ) {
193 initialize(x);
194 Real myval(0), ptval(0), val(0), one(1), two(2), error(two*tol + one);
195 std::vector<Real> ptvals;
196 while ( error > tol ) {
197 ValueSampler_->refine();
198 for ( int i = ValueSampler_->start(); i < ValueSampler_->numMySamples(); ++i ) {
199 getValue(ptval,x,ValueSampler_->getMyPoint(i),tol);
200 myval += ValueSampler_->getMyWeight(i)*ptval;
201 ptvals.push_back(ptval);
202 }
203 error = ValueSampler_->computeError(ptvals);
204 ptvals.clear();
205 }
206 ValueSampler_->sumAll(&myval,&val,1);
207 value_ += val;
208 ValueSampler_->setSamples();
209 tol = error;
210 return value_;
211 }
212
213 void gradient( Vector<Real> &g, const Vector<Real> &x, Real &tol ) {
214 initialize(x);
215 g.zero(); pointDual_->zero(); sumDual_->zero();
216 std::vector<Ptr<Vector<Real>>> ptgs;
217 Real one(1), two(2), error(two*tol + one);
218 while ( error > tol ) {
219 GradientSampler_->refine();
220 for ( int i = GradientSampler_->start(); i < GradientSampler_->numMySamples(); ++i ) {
221 getGradient(*pointDual_,x,GradientSampler_->getMyPoint(i),tol);
222 sumDual_->axpy(GradientSampler_->getMyWeight(i),*pointDual_);
223 ptgs.push_back(pointDual_->clone());
224 (ptgs.back())->set(*pointDual_);
225 }
226 error = GradientSampler_->computeError(ptgs,x);
227//if (GradientSampler_->batchID()==0) {
228// std::cout << "IN GRADIENT: ERROR = " << error << " TOL = " << tol << std::endl;
229//}
230 ptgs.clear();
231 }
232 GradientSampler_->sumAll(*sumDual_,g);
233 gradient_->plus(g);
234 g.set(*(gradient_));
235 GradientSampler_->setSamples();
236 tol = error;
237 }
238
239 void hessVec( Vector<Real> &hv, const Vector<Real> &v,
240 const Vector<Real> &x, Real &tol ) {
241 initialize(x);
242 hv.zero(); pointDual_->zero(); sumDual_->zero();
243 for ( int i = 0; i < HessianSampler_->numMySamples(); ++i ) {
244 getHessVec(*pointDual_,v,x,HessianSampler_->getMyPoint(i),tol);
245 sumDual_->axpy(HessianSampler_->getMyWeight(i),*pointDual_);
246 }
247 HessianSampler_->sumAll(*sumDual_,hv);
248 }
249
250 void precond( Vector<Real> &Pv, const Vector<Real> &v,
251 const Vector<Real> &x, Real &tol ) {
252 Pv.set(v.dual());
253 }
254};
255
256}
257
258#endif
Provides the interface to evaluate objective functions.
RiskNeutralObjective(const Ptr< Objective< Real > > &pObj, const Ptr< SampleGenerator< Real > > &sampler, const bool storage=true)
void hessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply Hessian approximation to vector.
Ptr< Objective< Real > > ParametrizedObjective_
Ptr< SampleGenerator< Real > > HessianSampler_
RiskNeutralObjective(const Ptr< Objective< Real > > &pObj, const Ptr< SampleGenerator< Real > > &vsampler, const Ptr< SampleGenerator< Real > > &gsampler, const bool storage=true)
void getHessVec(Vector< Real > &hv, const Vector< Real > &v, const Vector< Real > &x, const std::vector< Real > &param, Real &tol)
void update(const Vector< Real > &x, UpdateType type, int iter=-1)
Update objective function.
Ptr< SampleGenerator< Real > > GradientSampler_
void getValue(Real &val, const Vector< Real > &x, const std::vector< Real > &param, Real &tol)
void precond(Vector< Real > &Pv, const Vector< Real > &v, const Vector< Real > &x, Real &tol)
Apply preconditioner to vector.
RiskNeutralObjective(const Ptr< Objective< Real > > &pObj, const Ptr< SampleGenerator< Real > > &vsampler, const Ptr< SampleGenerator< Real > > &gsampler, const Ptr< SampleGenerator< Real > > &hsampler, const bool storage=true)
Real value(const Vector< Real > &x, Real &tol)
Compute value.
Ptr< SampleGenerator< Real > > ValueSampler_
Ptr< VectorController< Real > > gradient_storage_
Ptr< ScalarController< Real > > value_storage_
void initialize(const Vector< Real > &x)
void gradient(Vector< Real > &g, const Vector< Real > &x, Real &tol)
Compute gradient.
void update(const Vector< Real > &x, bool flag=true, int iter=-1)
Update objective function.
void getGradient(Vector< Real > &g, const Vector< Real > &x, const std::vector< Real > &param, Real &tol)
Defines the linear algebra or vector space interface.
virtual void set(const Vector &x)
Set where .
virtual const Vector & dual() const
Return dual representation of , for example, the result of applying a Riesz map, or change of basis,...
virtual void zero()
Set to zero vector.