ROL
ROL_MeanDeviation.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_MEANDEVIATION_HPP
45 #define ROL_MEANDEVIATION_HPP
46 
47 #include "ROL_RiskMeasure.hpp"
48 #include "ROL_PositiveFunction.hpp"
49 
50 namespace ROL {
51 
52 template<class Real>
53 class MeanDeviation : 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  MeanDeviation( 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  MeanDeviation( 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  Real diff = 0.0, pf0 = 0.0, dev = 0.0;
158  std::vector<Real> devp(this->order_.size(),0.0);
159  std::vector<Real> devs(this->order_.size(),0.0);
160  for ( unsigned i = 0; i < this->weights_.size(); i++ ) {
161  diff = this->value_storage_[i]-ev;
162  pf0 = this->positiveFunction_->evaluate(diff,0);
163  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
164  devp[p] += std::pow(pf0,this->order_[p]) * this->weights_[i];
165  }
166  }
167  sampler.sumAll(&devp[0],&devs[0],devp.size());
168  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
169  dev += (this->coeff_[p])*std::pow(devs[p],1.0/this->order_[p]);
170  }
171  // Return mean plus deviation
172  return ev + dev;
173  }
174 
176  g.zero();
177  // Compute expected value
178  Real val = RiskMeasure<Real>::val_;
179  Real ev = 0.0;
180  sampler.sumAll(&val,&ev,1);
181  // Compute deviation
182  Real diff = 0.0, pf0 = 0.0, pf1 = 0.0, c = 0.0;
183  std::vector<Real> dev0(this->order_.size(),0.0);
184  std::vector<Real> des0(this->order_.size(),0.0);
185  std::vector<Real> dev1(this->order_.size(),0.0);
186  std::vector<Real> des1(this->order_.size(),0.0);
187  for ( unsigned i = 0; i < this->weights_.size(); i++ ) {
188  diff = this->value_storage_[i]-ev;
189  pf0 = this->positiveFunction_->evaluate(diff,0);
190  pf1 = this->positiveFunction_->evaluate(diff,1);
191  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
192  dev0[p] += (this->weights_[i]) * std::pow(pf0,this->order_[p]);
193  dev1[p] += (this->weights_[i]) * std::pow(pf0,this->order_[p]-1.0) * pf1;
194  }
195  }
196  sampler.sumAll(&dev0[0],&des0[0],dev0.size());
197  sampler.sumAll(&dev1[0],&des1[0],dev1.size());
198  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
199  dev0[p] = std::pow(des0[p],1.0-1.0/this->order_[p]);
200  }
201  // Compute derivative
202  Teuchos::RCP<Vector<Real> > gtmp = g.clone();
203  gtmp->zero();
204  for ( unsigned i = 0; i < this->weights_.size(); i++ ) {
205  c = 0.0;
206  diff = this->value_storage_[i]-ev;
207  pf0 = this->positiveFunction_->evaluate(diff,0);
208  pf1 = this->positiveFunction_->evaluate(diff,1);
209  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
210  if ( dev0[p] > 0.0 ) {
211  c += (this->coeff_[p])/dev0[p] * (std::pow(pf0,this->order_[p]-1.0)*pf1 - des1[p]);
212  }
213  }
214  gtmp->axpy((this->weights_[i])*c,*(this->gradient_storage_[i]));
215  }
216  gtmp->axpy(1.0,*(RiskMeasure<Real>::g_));
217  sampler.sumAll(*gtmp,g);
218  }
219 
221  hv.zero();
222  // Compute expected value
223  Real val = RiskMeasure<Real>::val_;
224  Real ev = 0.0;
225  sampler.sumAll(&val,&ev,1);
226  Real gv = RiskMeasure<Real>::gv_;
227  Real egv = 0.0;
228  sampler.sumAll(&gv,&egv,1);
229  // Compute deviation
230  Real diff = 0.0, pf0 = 0.0, pf1 = 0.0, pf2 = 0.0;
231  Real cg = 0.0, ch = 0.0, diff1 = 0.0, diff2 = 0.0, diff3 = 0.0;
232  std::vector<Real> dev0(this->order_.size(),0.0);
233  std::vector<Real> dev1(this->order_.size(),0.0);
234  std::vector<Real> dev2(this->order_.size(),0.0);
235  std::vector<Real> dev3(this->order_.size(),0.0);
236  std::vector<Real> des0(this->order_.size(),0.0);
237  std::vector<Real> des1(this->order_.size(),0.0);
238  std::vector<Real> des2(this->order_.size(),0.0);
239  std::vector<Real> des3(this->order_.size(),0.0);
240  std::vector<Real> devp(this->order_.size(),0.0);
241  std::vector<Real> gvp1(this->order_.size(),0.0);
242  std::vector<Real> gvp2(this->order_.size(),0.0);
243  std::vector<Real> gvp3(this->order_.size(),0.0);
244  std::vector<Real> gvs1(this->order_.size(),0.0);
245  std::vector<Real> gvs2(this->order_.size(),0.0);
246  std::vector<Real> gvs3(this->order_.size(),0.0);
247  for ( unsigned i = 0; i < this->weights_.size(); i++ ) {
248  diff = this->value_storage_[i]-ev;
249  pf0 = (this->positiveFunction_)->evaluate(diff,0);
250  pf1 = (this->positiveFunction_)->evaluate(diff,1);
251  pf2 = (this->positiveFunction_)->evaluate(diff,2);
252  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
253  dev0[p] += (this->weights_[i]) * std::pow(pf0,this->order_[p]);
254  dev1[p] += (this->weights_[i]) * std::pow(pf0,this->order_[p]-1.0) * pf1;
255  dev2[p] += (this->weights_[i]) * std::pow(pf0,this->order_[p]-2.0) * pf1 * pf1;
256  dev3[p] += (this->weights_[i]) * std::pow(pf0,this->order_[p]-1.0) * pf2;
257  }
258  }
259  sampler.sumAll(&dev0[0],&des0[0],dev0.size());
260  sampler.sumAll(&dev1[0],&des1[0],dev1.size());
261  sampler.sumAll(&dev2[0],&des2[0],dev2.size());
262  sampler.sumAll(&dev3[0],&des3[0],dev3.size());
263  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
264  devp[p] = std::pow(des0[p],2.0-1.0/this->order_[p]);
265  dev0[p] = std::pow(des0[p],1.0-1.0/this->order_[p]);
266  }
267  for ( unsigned i = 0; i < this->value_storage_.size(); i++ ) {
268  diff = this->value_storage_[i]-ev;
269  pf0 = (this->positiveFunction_)->evaluate(diff,0);
270  pf1 = (this->positiveFunction_)->evaluate(diff,1);
271  pf2 = (this->positiveFunction_)->evaluate(diff,2);
272  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
273  gvp1[p] += (this->weights_[i]) * (std::pow(pf0,this->order_[p]-1.0)*pf1-des1[p]) *
274  (this->gradvec_storage_[i] - egv);
275  gvp2[p] += (this->weights_[i]) * (std::pow(pf0,this->order_[p]-2.0)*pf1*pf1-des2[p]) *
276  (this->gradvec_storage_[i] - egv);
277  gvp3[p] += (this->weights_[i]) * (std::pow(pf0,this->order_[p]-1.0)*pf2-des3[p]) *
278  (this->gradvec_storage_[i] - egv);
279  }
280  }
281  sampler.sumAll(&gvp1[0],&gvs1[0],gvp1.size());
282  sampler.sumAll(&gvp2[0],&gvs2[0],gvp2.size());
283  sampler.sumAll(&gvp3[0],&gvs3[0],gvp3.size());
284  // Compute derivative
285  Teuchos::RCP<Vector<Real> > htmp = hv.clone();
286  htmp->zero();
287  for ( unsigned i = 0; i < this->weights_.size(); i++ ) {
288  cg = 1.0;
289  ch = 0.0;
290  diff = this->value_storage_[i]-ev;
291  pf0 = (this->positiveFunction_)->evaluate(diff,0);
292  pf1 = (this->positiveFunction_)->evaluate(diff,1);
293  pf2 = (this->positiveFunction_)->evaluate(diff,2);
294  for ( unsigned p = 0; p < this->order_.size(); p++ ) {
295  if ( dev0[p] > 0.0 ) {
296  diff1 = std::pow(pf0,this->order_[p]-1.0)*pf1-des1[p];
297  diff2 = std::pow(pf0,this->order_[p]-2.0)*pf1*pf1*(this->gradvec_storage_[i]-egv)-gvs2[p];
298  diff3 = std::pow(pf0,this->order_[p]-1.0)*pf2*(this->gradvec_storage_[i]-egv)-gvs3[p];
299  cg += this->coeff_[p]*diff1/dev0[p];
300  ch += this->coeff_[p]*(((this->order_[p]-1.0)*diff2+diff3)/dev0[p] -
301  (this->order_[p]-1.0)*gvs1[p]*diff1/devp[p]);
302  }
303  }
304  htmp->axpy(this->weights_[i]*ch,*(this->gradient_storage_[i]));
305  htmp->axpy(this->weights_[i]*cg,*(this->hessvec_storage_[i]));
306  }
307  sampler.sumAll(*htmp,hv);
308  }
309 };
310 
311 }
312 
313 #endif
Real getValue(SampleGenerator< Real > &sampler)
MeanDeviation(Real order, Real coeff, Teuchos::RCP< PositiveFunction< Real > > &pf)
virtual Teuchos::RCP< Vector > clone() const =0
Clone to make a new (uninitialized) vector.
virtual void zero()
Set to zero vector.
Definition: ROL_Vector.hpp:155
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
void sumAll(Real *input, Real *output, int dim) const
void getGradient(Vector< Real > &g, 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 > order_
std::vector< Real > value_storage_
std::vector< Teuchos::RCP< Vector< Real > > > gradient_storage_
std::vector< Real > gradvec_storage_
void update(const Real val, const Real weight)
void update(const Real val, const Vector< Real > &g, const Real gv, const Vector< Real > &hv, const Real weight)
void update(const Real val, const Vector< Real > &g, const Real weight)
MeanDeviation(std::vector< Real > &order, std::vector< Real > &coeff, Teuchos::RCP< PositiveFunction< Real > > &pf)
Teuchos::RCP< PositiveFunction< Real > > positiveFunction_
void reset(Teuchos::RCP< Vector< Real > > &x0, const Vector< Real > &x)
std::vector< Real > weights_
std::vector< Real > coeff_
std::vector< Teuchos::RCP< Vector< Real > > > hessvec_storage_
void getHessVec(Vector< Real > &hv, SampleGenerator< Real > &sampler)