ROL
ROL_MeanVariance.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_MEANVARIANCE_HPP
45 #define ROL_MEANVARIANCE_HPP
46 
47 #include "ROL_RiskMeasure.hpp"
48 #include "ROL_PositiveFunction.hpp"
49 
50 namespace ROL {
51 
52 template<class Real>
53 class MeanVariance : public RiskMeasure<Real> {
54 private:
55  std::vector<Real> order_;
56  std::vector<Real> coeff_;
57 
58  std::vector<Real> weights_;
59  std::vector<Real> value_storage_;
60  std::vector<Teuchos::RCP<Vector<Real> > > gradient_storage_;
61  std::vector<Teuchos::RCP<Vector<Real> > > hessvec_storage_;
62  std::vector<Real> gradvec_storage_;
63 
64  Teuchos::RCP<PositiveFunction<Real> > positiveFunction_;
65 
66 public:
67  MeanVariance( Real order, Real coeff,
68  Teuchos::RCP<PositiveFunction<Real> > &pf ) : positiveFunction_(pf) {
69  order_.clear();
70  order_.push_back((order < 2.0) ? 2.0 : order);
71  coeff_.clear();
72  coeff_.push_back((coeff < 0.0) ? 1.0 : coeff);
73  }
74  MeanVariance( std::vector<Real> &order, std::vector<Real> &coeff,
75  Teuchos::RCP<PositiveFunction<Real> > &pf ) : positiveFunction_(pf) {
76  order_.clear();
77  coeff_.clear();
78  if ( order.size() != coeff.size() ) {
79  coeff.resize(order.size(),1.0);
80  }
81  for ( unsigned i = 0; i < order.size(); i++ ) {
82  order_.push_back((order[i] < 2.0) ? 2.0 : order[i]);
83  coeff_.push_back((coeff[i] < 0.0) ? 1.0 : coeff[i]);
84  }
85  }
86 
87  void reset(Teuchos::RCP<Vector<Real> > &x0, const Vector<Real> &x) {
92  (this->value_storage_).clear();
93  (this->gradient_storage_).clear();
94  (this->gradvec_storage_).clear();
95  (this->hessvec_storage_).clear();
96  (this->weights_).clear();
97  x0 = Teuchos::rcp(&const_cast<Vector<Real> &>(x),false);
98  }
99 
100  void reset(Teuchos::RCP<Vector<Real> > &x0, const Vector<Real> &x,
101  Teuchos::RCP<Vector<Real> > &v0, const Vector<Real> &v) {
106  (this->value_storage_).clear();
107  (this->gradient_storage_).clear();
108  (this->gradvec_storage_).clear();
109  (this->hessvec_storage_).clear();
110  (this->weights_).clear();
111  x0 = Teuchos::rcp(&const_cast<Vector<Real> &>(x),false);
112  v0 = Teuchos::rcp(&const_cast<Vector<Real> &>(v),false);
113  }
114 
115  void update(const Real val, const Real weight) {
116  RiskMeasure<Real>::val_ += weight * val;
117  this->value_storage_.push_back(val);
118  this->weights_.push_back(weight);
119  }
120 
121  void update(const Real val, const Vector<Real> &g, const Real weight) {
122  RiskMeasure<Real>::val_ += weight * val;
123  RiskMeasure<Real>::g_->axpy(weight,g);
124  this->value_storage_.push_back(val);
125  this->gradient_storage_.push_back(g.clone());
126  typename std::vector<Teuchos::RCP<Vector<Real> > >::iterator it = (this->gradient_storage_).end();
127  it--;
128  (*it)->set(g);
129  this->weights_.push_back(weight);
130  }
131 
132  void update(const Real val, const Vector<Real> &g, const Real gv, const Vector<Real> &hv,
133  const Real weight) {
134  RiskMeasure<Real>::val_ += weight * val;
135  RiskMeasure<Real>::gv_ += weight * gv;
136  RiskMeasure<Real>::g_->axpy(weight,g);
137  RiskMeasure<Real>::hv_->axpy(weight,hv);
138  this->value_storage_.push_back(val);
139  this->gradient_storage_.push_back(g.clone());
140  typename std::vector<Teuchos::RCP<Vector<Real> > >::iterator it = (this->gradient_storage_).end();
141  it--;
142  (*it)->set(g);
143  this->gradvec_storage_.push_back(gv);
144  this->hessvec_storage_.push_back(hv.clone());
145  it = (this->hessvec_storage_).end();
146  it--;
147  (*it)->set(hv);
148  this->weights_.push_back(weight);
149  }
150 
152  // Compute expected value
153  Real val = RiskMeasure<Real>::val_;
154  Real ev = 0.0;
155  sampler.sumAll(&val,&ev,1);
156  // Compute deviation
157  val = 0.0;
158  Real diff = 0.0, pf0 = 0.0, var = 0.0;
159  for ( unsigned i = 0; i < this->weights_.size(); i++ ) {
160  diff = this->value_storage_[i]-ev;
161  pf0 = this->positiveFunction_->evaluate(diff,0);
162  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
163  val += this->coeff_[p] * std::pow(pf0,this->order_[p]) * this->weights_[i];
164  }
165  }
166  sampler.sumAll(&val,&var,1);
167  // Return mean plus deviation
168  return ev + var;
169  }
170 
172  g.zero();
173  // Compute expected value
174  Real val = RiskMeasure<Real>::val_;
175  Real ev = 0.0;
176  sampler.sumAll(&val,&ev,1);
177  sampler.sumAll(*(RiskMeasure<Real>::g_),g);
178  // Compute deviation
179  Teuchos::RCP<Vector<Real> > gs = g.clone(); gs->zero();
180  Teuchos::RCP<Vector<Real> > gtmp = g.clone(); gtmp->zero();
181  Real diff = 0.0, pf0 = 0.0, pf1 = 0.0, c = 0.0, ec = 0.0, ecs = 0.0;
182  for ( unsigned i = 0; i < this->weights_.size(); i++ ) {
183  c = 0.0;
184  diff = this->value_storage_[i]-ev;
185  pf0 = this->positiveFunction_->evaluate(diff,0);
186  pf1 = this->positiveFunction_->evaluate(diff,1);
187  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
188  c += this->coeff_[p]*this->order_[p]*std::pow(pf0,this->order_[p]-1.0)*pf1;
189  }
190  ec += this->weights_[i]*c;
191  gtmp->axpy(this->weights_[i]*c,*(this->gradient_storage_[i]));
192  }
193  sampler.sumAll(&ec,&ecs,1);
194  g.scale(1.0-ecs);
195  sampler.sumAll(*gtmp,*gs);
196  g.plus(*gs);
197  }
198 
200  hv.zero();
201  // Compute expected value
202  Real val = RiskMeasure<Real>::val_;
203  Real ev = 0.0;
204  sampler.sumAll(&val,&ev,1);
205  Real gv = RiskMeasure<Real>::gv_;
206  Real egv = 0.0;
207  sampler.sumAll(&gv,&egv,1);
208  Teuchos::RCP<Vector<Real> > g = hv.clone();
209  sampler.sumAll(*(RiskMeasure<Real>::g_),*g);
210  sampler.sumAll(*(RiskMeasure<Real>::hv_),hv);
211  // Compute deviation
212  Real diff = 0.0, pf0 = 0.0, pf1 = 0.0, pf2 = 0.0;
213  Real cg = 0.0, ecg = 0.0, ecgs = 0.0, ch = 0.0, ech = 0.0, echs = 0.0;
214  Teuchos::RCP<Vector<Real> > htmp = hv.clone(); htmp->zero();
215  Teuchos::RCP<Vector<Real> > hs = hv.clone(); hs->zero();
216  for ( unsigned i = 0; i < this->weights_.size(); i++ ) {
217  cg = 0.0;
218  ch = 0.0;
219  diff = this->value_storage_[i]-ev;
220  pf0 = this->positiveFunction_->evaluate(diff,0);
221  pf1 = this->positiveFunction_->evaluate(diff,1);
222  pf2 = this->positiveFunction_->evaluate(diff,2);
223  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
224  cg += this->coeff_[p]*this->order_[p]*(this->gradvec_storage_[i]-egv)*
225  ((this->order_[p]-1.0)*std::pow(pf0,this->order_[p]-2.0)*pf1*pf1+
226  std::pow(pf0,this->order_[p]-1.0)*pf2);
227  ch += this->coeff_[p]*this->order_[p]*std::pow(pf0,this->order_[p]-1.0)*pf1;
228  }
229  ecg += this->weights_[i]*cg;
230  ech += this->weights_[i]*ch;
231  htmp->axpy(this->weights_[i]*cg,*(this->gradient_storage_[i]));
232  htmp->axpy(this->weights_[i]*ch,*(this->hessvec_storage_[i]));
233  }
234  sampler.sumAll(&ech,&echs,1);
235  hv.scale(1.0-echs);
236  sampler.sumAll(&ecg,&ecgs,1);
237  hv.axpy(-ecgs,*g);
238  sampler.sumAll(*htmp,*hs);
239  hv.plus(*hs);
240  }
241 };
242 
243 }
244 
245 #endif
Real getValue(SampleGenerator< Real > &sampler)
std::vector< Teuchos::RCP< Vector< Real > > > gradient_storage_
virtual void scale(const Real alpha)=0
Compute where .
MeanVariance(Real order, Real coeff, Teuchos::RCP< PositiveFunction< Real > > &pf)
std::vector< Real > gradvec_storage_
virtual void plus(const Vector &x)=0
Compute , where .
Teuchos::RCP< PositiveFunction< Real > > positiveFunction_
virtual void axpy(const Real alpha, const Vector &x)
Compute where .
Definition: ROL_Vector.hpp:141
void reset(Teuchos::RCP< Vector< Real > > &x0, const Vector< Real > &x)
void update(const Real val, const Real weight)
std::vector< Real > order_
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
std::vector< Real > weights_
MeanVariance(std::vector< Real > &order, std::vector< Real > &coeff, Teuchos::RCP< PositiveFunction< Real > > &pf)
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:155
std::vector< Real > coeff_
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
void update(const Real val, const Vector< Real > &g, const Real gv, const Vector< Real > &hv, const Real weight)
void sumAll(Real *input, Real *output, int dim) const
void update(const Real val, const Vector< Real > &g, const Real weight)
void getGradient(Vector< Real > &g, SampleGenerator< Real > &sampler)
std::vector< Teuchos::RCP< Vector< Real > > > hessvec_storage_
void getHessVec(Vector< Real > &hv, SampleGenerator< Real > &sampler)
void reset(Teuchos::RCP< Vector< Real > > &x0, const Vector< Real > &x, Teuchos::RCP< Vector< Real > > &v0, const Vector< Real > &v)
std::vector< Real > value_storage_