ROL
ROL_MonteCarloGenerator.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_MONTECARLOGENERATOR_HPP
45 #define ROL_MONTECARLOGENERATOR_HPP
46 
47 #include "ROL_SampleGenerator.hpp"
48 #include "ROL_Distribution.hpp"
49 
50 namespace ROL {
51 
52 template<class Real>
53 class MonteCarloGenerator : public SampleGenerator<Real> {
54 private:
55  int nSamp_;
57  bool use_SA_;
58  bool adaptive_;
60  std::vector<std::vector<Real> > data_;
61 
62  Real sum_val_;
63  Real sum_val2_;
64  Real sum_ng_;
65  Real sum_ng2_;
66 
67  bool useDist_;
68  std::vector<Teuchos::RCP<ROL::Distribution<Real> > > dist_;
69 
70  Real ierf(Real input) {
71  std::vector<Real> coeff;
72  Real c = 1.0;
73  Real tmp = c * (std::sqrt(M_PI)/2.0 * input);
74  Real val = tmp;
75  coeff.push_back(c);
76  int cnt = 1;
77  while (std::abs(tmp) > 1.e-4*std::abs(val)) {
78  c = 0.0;
79  for ( unsigned i = 0; i < coeff.size(); i++ ) {
80  c += coeff[i]*coeff[coeff.size()-1-i]/((i+1)*(2*i+1));
81  }
82  tmp = c/(2.0*(Real)cnt+1.0) * std::pow(std::sqrt(M_PI)/2.0 * input,2.0*(Real)cnt+1.0);
83  val += tmp;
84  coeff.push_back(c);
85  cnt++;
86  }
87  return val;
88  }
89 
90  void sample(void) {
91  // Get process rank and number of processes
92  int rank = SampleGenerator<Real>::batchID();
94  // Separate samples across processes
95  int frac = this->nSamp_ / nProc;
96  int rem = this->nSamp_ % nProc;
97  unsigned N = frac;
98  if ( rank < rem ) {
99  N++;
100  }
101  // Generate samples
102  std::vector<std::vector<Real> > pts;
103  std::vector<Real> p;
104  //srand((rank+1)*(rank+1)*time(NULL));
105  for ( unsigned i = 0; i < N; i++ ) {
106  srand(123456*(rank*N + (i+1)));
107  if ( !useDist_ ) {
108  p.resize(this->data_.size(),0.0);
109  for ( unsigned j = 0; j < this->data_.size(); j++ ) {
110  if ( this->use_normal_ ) {
111  p[j] = std::sqrt(2.0*(this->data_[j])[1])*this->ierf(2.0*((Real)rand())/((Real)RAND_MAX)-1.0) +
112  (this->data_[j])[0];
113  }
114  else {
115  p[j] = ((this->data_[j])[1]-(this->data_[j])[0])*((Real)rand())/((Real)RAND_MAX)+(this->data_[j])[0];
116  }
117  }
118  }
119  else {
120  p.resize(dist_.size(),0.0);
121  for ( unsigned j = 0; j < dist_.size(); j++ ) {
122  p[j] = (dist_[j])->invcdf((Real)rand()/(Real)RAND_MAX);
123  while (std::abs(p[j]) > 0.1*ROL::ROL_OVERFLOW) {
124  p[j] = (dist_[j])->invcdf((Real)rand()/(Real)RAND_MAX);
125  }
126  }
127  }
128  pts.push_back(p);
129  }
130  std::vector<Real> wts(N,1.0/((Real)this->nSamp_));
133  }
134 
135  std::vector<std::vector<Real> > sample(int nSamp, bool store = true) {
136  // Get process rank and number of processes
137  int rank = SampleGenerator<Real>::batchID();
138  int nProc = SampleGenerator<Real>::numBatches();
139  // Separate samples across processes
140  int frac = nSamp / nProc;
141  int rem = nSamp % nProc;
142  unsigned N = frac;
143  if ( rank < rem ) {
144  N++;
145  }
146  // Generate samples
147  std::vector<std::vector<Real> > pts;
148  std::vector<Real> p;
149  //srand((rank+1)*(rank+1)*time(NULL));
150  for ( unsigned i = 0; i < N; i++ ) {
151  srand(123456*(rank*N + (i+1)));
152  if ( !useDist_ ) {
153  p.resize(this->data_.size(),0.0);
154  for ( unsigned j = 0; j < this->data_.size(); j++ ) {
155  if ( this->use_normal_ ) {
156  p[j] = std::sqrt(2.0*(this->data_[j])[1])*this->ierf(2.0*((Real)rand())/((Real)RAND_MAX)-1.0) +
157  (this->data_[j])[0];
158  }
159  else {
160  p[j] = ((this->data_[j])[1]-(this->data_[j])[0])*((Real)rand())/((Real)RAND_MAX)+(this->data_[j])[0];
161  }
162  }
163  }
164  else {
165  p.resize(dist_.size(),0.0);
166  for ( unsigned j = 0; j < dist_.size(); j++ ) {
167  p[j] = (dist_[j])->invcdf((Real)rand()/(Real)RAND_MAX);
168  while (std::abs(p[j]) > 0.1*ROL::ROL_OVERFLOW) {
169  p[j] = (dist_[j])->invcdf((Real)rand()/(Real)RAND_MAX);
170  }
171  }
172  }
173  pts.push_back(p);
174  }
175  if ( store ) {
176  std::vector<Real> wts(N,1.0/((Real)nSamp));
179  }
180  return pts;
181  }
182 
183 public:
184  MonteCarloGenerator( int nSamp, std::vector<Teuchos::RCP<Distribution<Real> > > &dist,
185  Teuchos::RCP<BatchManager<Real> > &bman,
186  bool use_SA = false, bool adaptive = false, int numNewSamps = 0 )
187  : SampleGenerator<Real>(bman), nSamp_(nSamp), use_normal_(false), use_SA_(use_SA), adaptive_(adaptive),
188  numNewSamps_(numNewSamps), sum_val_(0.0), sum_val2_(0.0), sum_ng_(0.0), sum_ng2_(0.0),
189  useDist_(true), dist_(dist) {
190  sample();
191  }
192 
193  MonteCarloGenerator( int nSamp, std::vector<std::vector<Real> > &bounds,
194  Teuchos::RCP<BatchManager<Real> > &bman,
195  bool use_SA = false, bool adaptive = false, int numNewSamps = 0 )
196  : SampleGenerator<Real>(bman), nSamp_(nSamp), use_normal_(false), use_SA_(use_SA), adaptive_(adaptive),
197  numNewSamps_(numNewSamps), sum_val_(0.0), sum_val2_(0.0), sum_ng_(0.0), sum_ng2_(0.0),
198  useDist_(false) {
199  unsigned dim = bounds.size();
200  data_.clear();
201  Real tmp = 0.0;
202  for ( unsigned j = 0; j < dim; j++ ) {
203  if ( (bounds[j])[0] > (bounds[j])[1] ) {
204  tmp = (bounds[j])[0];
205  (bounds[j])[0] = (bounds[j])[1];
206  (bounds[j])[1] = tmp;
207  data_.push_back(bounds[j]);
208  }
209  data_.push_back(bounds[j]);
210  }
211  sample();
212  }
213 
214  MonteCarloGenerator( int nSamp, std::vector<Real> mean, std::vector<Real> std,
215  Teuchos::RCP<BatchManager<Real> > &bman,
216  bool use_SA = false, bool adaptive = false, int numNewSamps = 0 )
217  : SampleGenerator<Real>(bman), nSamp_(nSamp), use_normal_(true), use_SA_(use_SA), adaptive_(adaptive),
218  numNewSamps_(numNewSamps), sum_val_(0.0), sum_val2_(0.0), sum_ng_(0.0), sum_ng2_(0.0),
219  useDist_(false) {
220  unsigned dim = mean.size();
221  data_.clear();
222  std::vector<Real> tmp(2,0.0);
223  for ( unsigned j = 0; j < dim; j++ ) {
224  tmp[0] = mean[j];
225  tmp[1] = std[j];
226  data_.push_back(tmp);
227  }
228  sample();
229  }
230 
231  void update( const Vector<Real> &x ) {
233  this->sum_val_ = 0.0;
234  this->sum_val2_ = 0.0;
235  this->sum_ng_ = 0.0;
236  this->sum_ng_ = 0.0;
237  if ( this->use_SA_ ) {
238  this->sample();
239  }
240  }
241 
242  Real computeError( std::vector<Real> &vals ) {
243  if ( this->adaptive_ && !(this->use_SA_) ) {
244  // Compute unbiased sample variance
245  int cnt = 0;
246  for ( int i = SampleGenerator<Real>::start(); i < SampleGenerator<Real>::numMySamples(); i++ ) {
247  this->sum_val_ += vals[cnt];
248  this->sum_val2_ += vals[cnt]*vals[cnt];
249  cnt++;
250  }
251  Real mymean = this->sum_val_ / this->nSamp_;
252  Real mean = 0.0;
253  SampleGenerator<Real>::sumAll(&mymean,&mean,1);
254 
255  Real myvar = (this->sum_val2_ - mean*mean)/(this->nSamp_-1.0);
256  Real var = 0.0;
257  SampleGenerator<Real>::sumAll(&myvar,&var,1);
258  // Return Monte Carlo error
259  vals.clear();
260  return std::sqrt(var/(this->nSamp_))*1.e-8;
261  }
262  else {
263  vals.clear();
264  return 0.0;
265  }
266  }
267 
268  Real computeError( std::vector<Teuchos::RCP<Vector<Real> > > &vals, const Vector<Real> &x ) {
269  if ( this->adaptive_ && !(this->use_SA_) ) {
270  // Compute unbiased sample variance
271  int cnt = 0;
272  Real ng = 0.0;
273  for ( int i = SampleGenerator<Real>::start(); i < SampleGenerator<Real>::numMySamples(); i++ ) {
274  ng = (vals[cnt])->norm();
275  this->sum_ng_ += ng;
276  this->sum_ng2_ += ng*ng;
277  cnt++;
278  }
279  Real mymean = this->sum_ng_ / this->nSamp_;
280  Real mean = 0.0;
281  SampleGenerator<Real>::sumAll(&mymean,&mean,1);
282 
283  Real myvar = (this->sum_ng2_ - mean*mean)/(this->nSamp_-1.0);
284  Real var = 0.0;
285  SampleGenerator<Real>::sumAll(&myvar,&var,1);
286  // Return Monte Carlo error
287  vals.clear();
288  return std::sqrt(var/(this->nSamp_))*1.e-4;
289  }
290  else {
291  vals.clear();
292  return 0.0;
293  }
294  }
295 
296  void refine(void) {
297  if ( this->adaptive_ && !(this->use_SA_) ) {
298  std::vector<std::vector<Real> > pts;
299  std::vector<Real> pt(this->data_.size(),0.0);
300  for ( int i = 0; i < SampleGenerator<Real>::numMySamples(); i++ ) {
302  pts.push_back(pt);
303  }
304  std::vector<std::vector<Real> > pts_new = this->sample(this->numNewSamps_,false);
305  pts.insert(pts.end(),pts_new.begin(),pts_new.end());
306  this->nSamp_ += this->numNewSamps_;
307  std::vector<Real> wts(pts.size(),1.0/((Real)this->nSamp_));
311  }
312  }
313 
314 };
315 
316 }
317 
318 #endif
virtual void update(const Vector< Real > &x)
virtual std::vector< Real > getMyPoint(const int i) const
Real computeError(std::vector< Teuchos::RCP< Vector< Real > > > &vals, const Vector< Real > &x)
Defines the linear algebra or vector space interface.
Definition: ROL_Vector.hpp:72
void sumAll(Real *input, Real *output, int dim) const
Real computeError(std::vector< Real > &vals)
std::vector< std::vector< Real > > sample(int nSamp, bool store=true)
virtual void refine(void)
MonteCarloGenerator(int nSamp, std::vector< std::vector< Real > > &bounds, Teuchos::RCP< BatchManager< Real > > &bman, bool use_SA=false, bool adaptive=false, int numNewSamps=0)
std::vector< Teuchos::RCP< ROL::Distribution< Real > > > dist_
MonteCarloGenerator(int nSamp, std::vector< Teuchos::RCP< Distribution< Real > > > &dist, Teuchos::RCP< BatchManager< Real > > &bman, bool use_SA=false, bool adaptive=false, int numNewSamps=0)
MonteCarloGenerator(int nSamp, std::vector< Real > mean, std::vector< Real > std, Teuchos::RCP< BatchManager< Real > > &bman, bool use_SA=false, bool adaptive=false, int numNewSamps=0)
void setPoints(std::vector< std::vector< Real > > &p)
static const double ROL_OVERFLOW
Platform-dependent maximum double.
Definition: ROL_Types.hpp:123
void update(const Vector< Real > &x)
void setWeights(std::vector< Real > &w)
std::vector< std::vector< Real > > data_