All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
mcarlo_vegas.h
Go to the documentation of this file.
1 /*
2  -------------------------------------------------------------------
3 
4  Copyright (C) 2006-2014, Andrew W. Steiner
5 
6  This file is part of O2scl.
7 
8  O2scl is free software; you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version.
12 
13  O2scl is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with O2scl. If not, see <http://www.gnu.org/licenses/>.
20 
21  -------------------------------------------------------------------
22 */
23 /* monte/vegas.c
24  *
25  * Copyright (C) 1996, 1997, 1998, 1999, 2000 Michael Booth
26  *
27  * This program is free software; you can redistribute it and/or modify
28  * it under the terms of the GNU General Public License as published by
29  * the Free Software Foundation; either version 3 of the License, or (at
30  * your option) any later version.
31  *
32  * This program is distributed in the hope that it will be useful, but
33  * WITHOUT ANY WARRANTY; without even the implied warranty of
34  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
35  * General Public License for more details.
36  *
37  * You should have received a copy of the GNU General Public License
38  * along with this program; if not, write to the Free Software
39  * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
40  * 02110-1301, USA.
41  */
42 #ifndef O2SCL_MCARLO_VEGAS_H
43 #define O2SCL_MCARLO_VEGAS_H
44 
45 /** \file mcarlo_vegas.h
46  \brief File defining \ref o2scl::mcarlo_vegas
47 */
48 
49 #include <iostream>
50 
51 #include <gsl/gsl_math.h>
52 #include <gsl/gsl_monte.h>
53 #include <gsl/gsl_machine.h>
54 #include <gsl/gsl_monte_vegas.h>
55 
56 #include <o2scl/mcarlo.h>
57 
58 #ifndef DOXYGEN_NO_O2NS
59 namespace o2scl {
60 #endif
61 
62  /** \brief Multidimensional integration using Vegas Monte Carlo (GSL)
63 
64  The output options are a little different than the original GSL
65  routine. The default setting of \ref mcarlo::verbose is 0,
66  which turns off all output. A verbose value of 1 prints summary
67  information about the weighted average and final result, while a
68  value of 2 also displays the grid coordinates. A value of 3
69  prints information from the rebinning procedure for each
70  iteration.
71 
72  Some original documentation from GSL:
73 
74  \verbatim
75  The input coordinates are x[j], with upper and lower limits
76  xu[j] and xl[j]. The integration length in the j-th direction is
77  delx[j]. Each coordinate x[j] is rescaled to a variable y[j] in
78  the range 0 to 1. The range is divided into bins with boundaries
79  xi[i][j], where i=0 corresponds to y=0 and i=bins to y=1. The
80  grid is refined (ie, bins are adjusted) using d[i][j] which is
81  some variation on the squared sum. A third parameter used in
82  defining the real coordinate using random numbers is called z.
83  It ranges from 0 to bins. Its integer part gives the lower index
84  of the bin into which a call is to be placed, and the remainder
85  gives the location inside the bin.
86 
87  When stratified sampling is used the bins are grouped into
88  boxes, and the algorithm allocates an equal number of function
89  calls to each box.
90 
91  The variable alpha controls how "stiff" the rebinning algorithm
92  is. alpha = 0 means never change the grid. Alpha is typically
93  set between 1 and 2.
94  \endverbatim
95 
96  \todo Mode = importance only doesn't give the same answer
97  as GSL yet.
98 
99  \future Prettify the verbose output
100 
101  \future Allow the user to get information about the how
102  the sampling was done, possibly by converting the bins
103  and boxes into a structure or class.
104 
105  \future Allow the user to change the maximum number of bins.
106 
107  Based on \ref Lepage78 . The current version of the algorithm
108  was described in the Cornell preprint CLNS-80/447 of March,
109  1980. The GSL code follows most closely the C version by D. R.
110  Yennie, coded in 1984.
111  */
112  template<class func_t=multi_funct11,
114  class rng_t=int, class rng_dist_t=rng_gsl>
115  //class rng_t=std::mt19937,
116  // class rng_dist_t=std::uniform_real_distribution<double> >
117  class mcarlo_vegas : public mcarlo<func_t,vec_t,rng_t,rng_dist_t> {
118 
119  public:
120 
124 
125  /// \name Integration mode (default is mode_importance)
126  //@{
127  int mode;
128  static const int mode_importance=1;
129  static const int mode_importance_only=0;
130  static const int mode_stratified=-1;
131  //@}
132 
133  /// Result from last iteration
134  double result;
135 
136  /// Uncertainty from last iteration
137  double sigma;
138 
139  /** \brief The stiffness of the rebinning algorithm (default 1.5)
140 
141  This usual range is between 1 and 2.
142  */
143  double alpha;
144 
145  /// Set the number of iterations (default 5)
146  unsigned int iterations;
147 
148  /** \brief The chi-squared per degree of freedom for the weighted
149  estimate of the integral
150 
151  After an integration, this should be close to 1. If it is not,
152  then this indicates that the values of the integral from
153  different iterations are inconsistent, and the error may be
154  underestimated. Further iterations of the algorithm may enable
155  one to obtain more reliable results.
156  */
157  double chisq;
158 
159  /// The output stream to send output information (default \c std::cout)
160  std::ostream *outs;
161 
162 #ifndef DOXYGEN_INTERNAL
163 
164  protected:
165 
166  /// Maximum number of bins
167  static const size_t bins_max=50;
168 
169  /// Number of dimensions
170  size_t dim;
171  /// Number of bins
172  unsigned int bins;
173  /// Number of boxes
174  unsigned int boxes;
175  /// Boundaries for each bin
177  /// Storage for grid refinement
179  /// The iteration length in each direction
181  /// The weight for each bin
183  /// The volume of the current bin
184  double vol;
185 
186  /// The bins for each direction
188  /// The boxes for each direction
190 
191  /// Distribution
193 
194  /** \name Scratch variables preserved between calls to
195  vegas_minteg_err()
196  */
197  //@{
198  double jac;
199  double wtd_int_sum;
200  double sum_wgts;
201  double chi_sum;
202  //@}
203 
204  /// The starting iteration number
205  unsigned int it_start;
206  /// The total number of iterations
207  unsigned int it_num;
208  /// Number of samples for computing chi squared
209  unsigned int samples;
210  /// Number of function calls per box
211  unsigned int calls_per_box;
212 
213  /// Initialize box coordinates
214  virtual void init_box_coord(ubvector_int &boxt) {
215  size_t i;
216  for (i=0;i<dim;i++) {
217  boxt[i]=0;
218  }
219  return;
220  }
221 
222  /** \brief Change box coordinates
223 
224  Steps through the box coordinates, e.g.
225  \verbatim
226  {0,0}, {0,1}, {0,2}, {0,3}, {1,0}, {1,1}, {1,2}, ...
227  \endverbatim
228 
229  This is among the member functions that is not virtual
230  because it is part of the innermost loop.
231  */
233  int j=dim-1;
234  int ng=boxes;
235 
236  while (j >= 0) {
237 
238  boxt[j]=(boxt[j]+1) % ng;
239  if (boxt[j] != 0) {
240  return 1;
241  }
242  j--;
243 
244  }
245 
246  return 0;
247  }
248 
249  /// Initialize grid
250  virtual void init_grid(const vec_t &xl, const vec_t &xu, size_t ldim) {
251  size_t j;
252  vol=1.0;
253  bins=1;
254 
255  for (j=0;j<ldim;j++) {
256  double dx=xu[j]-xl[j];
257  delx[j]=dx;
258  vol *= dx;
259 
260  xi[j]=0.0;
261  xi[dim+j]=1.0;
262  }
263 
264  return;
265  }
266 
267  /// Reset grid values
268  virtual void reset_grid_values() {
269  size_t i, j;
270 
271  for (i=0;i<bins;i++) {
272  for (j=0;j<dim;j++) {
273  d[i*dim+j]=0.0;
274  }
275  }
276  return;
277  }
278 
279  /** \brief Add the most recently generated result to the distribution
280 
281  This is among the member functions that is not virtual
282  because it is part of the innermost loop.
283  */
284  void accumulate_distribution(ubvector_int &lbin, double y) {
285  size_t j;
286 
287  for (j=0;j<dim;j++) {
288  int i=lbin[j];
289  d[i*dim+j]+=y;
290  }
291  return;
292  }
293 
294  /** \brief Generate a random position in a given box
295 
296  Use the random number generator to return a random position x
297  in a given box. The value of bin gives the bin location of the
298  random position (there may be several bins within a given box)
299 
300  This is among the member functions that is not virtual
301  because it is part of the innermost loop.
302  */
303  void random_point(vec_t &lx, ubvector_int &lbin, double &bin_vol,
304  const ubvector_int &lbox, const vec_t &xl,
305  const vec_t &xu) {
306 
307  double lvol=1.0;
308 
309  size_t j;
310 
311  for (j=0;j<dim;++j) {
312 
313  // The equivalent of gsl_rng_uniform_pos()
314  double rdn;
315  do {
316  rdn=this->rng_dist(this->rng);
317  } while (rdn==0);
318 
319  /* lbox[j] + ran gives the position in the box units, while z
320  is the position in bin units. */
321  double z=((lbox[j]+rdn)/boxes)*bins;
322 
323  int k=(int)z;
324 
325  double y, bin_width;
326 
327  lbin[j]=k;
328 
329  if (k == 0) {
330  bin_width=xi[dim+j];
331  y=z*bin_width;
332  } else {
333  bin_width=xi[(k+1)*dim+j]-xi[k*dim+j];
334  y=xi[k*dim+j]+(z-k)*bin_width;
335  }
336 
337  lx[j]=xl[j]+y*delx[j];
338 
339  lvol *= bin_width;
340  }
341 
342  bin_vol=lvol;
343 
344  return;
345  }
346 
347  /// Resize the grid
348  virtual void resize_grid(unsigned int lbins) {
349  size_t j, k;
350 
351  /* weight is ratio of bin sizes */
352 
353  double pts_per_bin=(double) bins/(double) lbins;
354 
355  for (j=0;j<dim;j++) {
356  double xold;
357  double xnew=0;
358  double dw=0;
359  int i=1;
360 
361  for (k=1;k <= bins;k++) {
362  dw+=1.0;
363  xold=xnew;
364  xnew=xi[k*dim+j];
365 
366  for (;dw > pts_per_bin;i++) {
367  dw -= pts_per_bin;
368  (xin[i])=xnew-(xnew-xold)*dw;
369  }
370  }
371 
372  for (k=1 ;k<lbins;k++) {
373  xi[k*dim+j]=xin[k];
374  }
375 
376  xi[lbins*dim+j]=1;
377  }
378 
379  bins=lbins;
380 
381  return;
382  }
383 
384  /// Refine the grid
385  virtual void refine_grid() {
386 
387  size_t i, j, k;
388 
389  for (j=0;j<dim;j++) {
390  double grid_tot_j, tot_weight;
391 
392  double oldg=d[j];
393  double newg=d[dim+j];
394 
395  d[j]=(oldg+newg)/2;
396  grid_tot_j=d[j];
397 
398  /* This implements gs[i][j]=(gs[i-1][j]+gs[i][j]+gs[i+1][j])/3 */
399 
400  for (i=1;i<bins-1;i++) {
401  double rc=oldg+newg;
402  oldg=newg;
403  newg=d[(i+1)*dim+j];
404  d[i*dim+j]=(rc+newg)/3;
405  grid_tot_j+=d[i*dim+j];
406  }
407  d[(bins-1)*dim+j]=(newg+oldg)/2;
408 
409  grid_tot_j+=d[(bins-1)*dim+j];
410 
411  tot_weight=0;
412 
413  for (i=0;i<bins;i++) {
414  weight[i]=0;
415 
416  if (d[i*dim+j] > 0) {
417  oldg=grid_tot_j/d[i*dim+j];
418  /* damped change */
419  weight[i]=pow (((oldg-1)/oldg/log (oldg)), alpha);
420  }
421 
422  tot_weight+=weight[i];
423  }
424 
425  {
426  double pts_per_bin=tot_weight/bins;
427 
428  double xold;
429  double xnew=0;
430  double dw=0;
431  i=1;
432 
433  for (k=0;k<bins;k++) {
434  dw+=weight[k];
435  xold=xnew;
436  xnew=xi[(k+1)*dim+j];
437 
438  for (;dw > pts_per_bin;i++) {
439  dw -= pts_per_bin;
440  (xin[i])=xnew-(xnew-xold)*dw/weight[k];
441  }
442  }
443 
444  for (k=1 ;k<bins ;k++) {
445  xi[k*dim+j]=xin[k];
446  }
447 
448  xi[bins*dim+j]=1;
449  }
450  }
451  return;
452  }
453 
454  /// Print limits of integration
455  virtual void print_lim(const vec_t &xl, const vec_t &xu,
456  unsigned long ldim) {
457 
458  unsigned long j;
459 
460  (*outs) << "The limits of integration are:\n" << std::endl;
461  for (j=0;j<ldim;++j) {
462  (*outs) << "xl[" << j << "]=" << xl[j] << " xu[" << j
463  << "]=" << xu[j] << std::endl;
464  }
465  (*outs) << std::endl;
466 
467  return;
468  }
469 
470  /// Print header
471  virtual void print_head(unsigned long num_dim, unsigned long calls,
472  unsigned int lit_num, unsigned int lbins,
473  unsigned int lboxes) {
474 
475  (*outs) << "num_dim=" << num_dim << ", calls=" << calls
476  << ", it_num=" << lit_num << ", max_it_num="
477  << iterations << ", verb=" << this->verbose << ", alph=" << alpha
478  << ",\n mode=" << mode << ", boxes=" << lboxes
479  << "\n\n single.......iteration "
480  << "accumulated......results\n"
481  << "iter integral sigma integral "
482  << " sigma chi-sq/it" << std::endl;
483 
484  return;
485  }
486 
487  /// Print results
488  virtual void print_res(unsigned int itr, double res,
489  double err, double cum_res, double cum_err,
490  double chi_sq) {
491  (*outs).precision(5);
492  (*outs) << itr << " ";
493  outs->setf(std::ios::showpos);
494  (*outs) << res << " ";
495  outs->unsetf(std::ios::showpos);
496  (*outs) << err << " ";
497  outs->setf(std::ios::showpos);
498  (*outs) << cum_res << " ";
499  outs->unsetf(std::ios::showpos);
500  (*outs) << cum_err << " " << chi_sq << std::endl;
501  (*outs).precision(6);
502  return;
503  }
504 
505  /// Print distribution
506  virtual void print_dist(unsigned long ldim) {
507  unsigned long i, j;
508 
509  if (this->verbose<2) {
510  return;
511  }
512 
513  for (j=0;j<ldim;++j) {
514  (*outs) << "\n Axis " << j << std::endl;
515  (*outs) << " x g" << std::endl;
516  outs->setf(std::ios::showpos);
517  for (i=0;i<bins;i++) {
518  (*outs) << "weight [ " << (xi[(i)*dim+(j)]) << " , "
519  << xi[(i+1)*dim+j] << " ] = ";
520  (*outs) << " " << (d[(i)*dim+(j)]) << std::endl;
521  }
522  outs->unsetf(std::ios::showpos);
523  (*outs) << std::endl;
524  }
525  (*outs) << std::endl;
526  return;
527  }
528 
529  /// Print grid
530  virtual void print_grid(unsigned long ldim) {
531 
532  if (this->verbose<2) {
533  return;
534  }
535 
536  unsigned long i, j;
537  for (j=0;j<ldim;++j) {
538  (*outs) << "\n Axis " << j << std::endl;
539  (*outs) << " x " << std::endl;
540  outs->setf(std::ios::showpos);
541  for (i=0;i <= bins;i++) {
542  (*outs) << (xi[(i)*dim+(j)]) << " ";
543  if (i % 5 == 4) {
544  (*outs) << std::endl;
545  }
546  }
547  (*outs) << std::endl;
548  outs->unsetf(std::ios::showpos);
549  }
550  (*outs) << std::endl;
551  return;
552  }
553 
554  /// Point for function evaluation
555  vec_t x;
556 
557 #endif
558 
559  public:
560 
561  mcarlo_vegas() {
562  this->verbose=0;
563  outs=&std::cout;
564  alpha=1.5;
565  iterations=5;
566  mode=mode_importance;
567  chisq=0;
568  bins=bins_max;
569  dim=0;
570  }
571 
572  /// Allocate memory
573  virtual int allocate(size_t ldim) {
574 
575  delx.resize(ldim);
576  d.resize(bins_max*ldim);
577  xi.resize((bins_max+1)*ldim);
578  xin.resize(bins_max+1);
579  weight.resize(bins_max);
580  box.resize(ldim);
581  bin.resize(ldim);
582  x.resize(ldim);
583 
584  dim=ldim;
585 
586  return 0;
587  }
588 
589  /** \brief Integrate function \c func from x=a to x=b.
590 
591  Original documentation from GSL:
592 
593  Normally, <tt>stage = 0</tt> which begins with a new uniform
594  grid and empty weighted average. Calling vegas with <tt>stage
595  = 1</tt> retains the grid from the previous run but discards
596  the weighted average, so that one can "tune" the grid using a
597  relatively small number of points and then do a large run with
598  <tt>stage = 1</tt> on the optimized grid. Setting <tt>stage =
599  2</tt> keeps the grid and the weighted average from the
600  previous run, but may increase (or decrease) the number of
601  histogram bins in the grid depending on the number of calls
602  available. Choosing <tt>stage = 3</tt> enters at the main
603  loop, so that nothing is changed, and is equivalent to
604  performing additional iterations in a previous call.
605 
606  \todo Should stage be passed by reference?
607  \todo There was an update between gsl-1.12 and 1.15 which
608  has not been implemented here yet.
609  */
610  virtual int vegas_minteg_err(int stage, func_t &func, size_t ndim,
611  const vec_t &xl, const vec_t &xu,
612  double &res, double &err) {
613 
614  size_t calls=this->n_points;
615 
616  double cum_int, cum_sig;
617  size_t i, k, it;
618 
619  for (i=0;i<dim;i++) {
620  if (xu[i] <= xl[i]) {
621  std::string serr="Upper limit, "+dtos(xu[i])+", must be greater "+
622  "than lower limit, "+dtos(xl[i])+", in mcarlo_vegas::"+
623  "vegas_minteg_err().";
624  O2SCL_ERR(serr.c_str(),exc_einval);
625  }
626 
627  if (xu[i]-xl[i] > GSL_DBL_MAX) {
628  O2SCL_ERR2("Range of integration is too large, please rescale ",
629  "in mcarlo_vegas::vegas_minteg_err().",exc_einval);
630  }
631  }
632 
633  if (stage == 0) {
634  init_grid(xl,xu,dim);
635  if (this->verbose>=1) {
636  print_lim(xl,xu,dim);
637  }
638  }
639 
640  if (stage<=1) {
641  wtd_int_sum=0;
642  sum_wgts=0;
643  chi_sum=0;
644  it_num=1;
645  samples=0;
646  chisq=0;
647  }
648 
649  if (stage <= 2) {
650 
651  unsigned int lbins=bins_max;
652  unsigned int lboxes=1;
653 
654  if (mode != mode_importance_only) {
655 
656  /* shooting for 2 calls/box */
657 
658  // The original GSL code was:
659  // boxes=floor (pow (calls/2.0, 1.0/dim));
660  // but floor returns double on my machine, so
661  // we explicitly typecast here
662 
663  lboxes=((unsigned int)(floor(pow(calls/2.0,1.0/dim))));
664  mode=mode_importance;
665 
666  if (2*lboxes >= bins_max) {
667  /* if bins/box < 2 */
668  int box_per_bin=GSL_MAX(lboxes/bins_max,1);
669 
670  if (lboxes/box_per_bin<bins_max) lbins=lboxes/box_per_bin;
671  else lbins=bins_max;
672  lboxes=box_per_bin*lbins;
673 
674  mode=mode_stratified;
675  }
676 
677  }
678 
679  //double tot_boxes=gsl_pow_int((double)boxes,dim);
680  double tot_boxes=pow((double)lboxes,(double)dim);
681  calls_per_box=((unsigned int)(GSL_MAX(calls/tot_boxes,2)));
682  calls=((size_t)( calls_per_box*tot_boxes));
683 
684  /* total volume of x-space/(avg num of calls/bin) */
685  jac=vol*pow((double) lbins, (double) dim)/calls;
686 
687  boxes=lboxes;
688 
689  /* If the number of bins changes from the previous invocation, bins
690  are expanded or contracted accordingly, while preserving bin
691  density */
692 
693  if (lbins!=bins) {
694  resize_grid(lbins);
695  if (this->verbose > 2) print_grid(dim);
696  }
697  if (this->verbose >= 1) {
698  print_head(dim,calls,it_num,bins,boxes);
699  }
700  }
701 
703 
704  cum_int=0.0;
705  cum_sig=0.0;
706 
707  for (it=0;it<iterations;it++) {
708 
709  double intgrl=0.0, intgrl_sq=0.0;
710  double tss=0.0;
711  double wgt, var, sig;
712  size_t lcalls_per_box=calls_per_box;
713  double jacbin=jac;
714 
715  it_num=it_start+it;
716 
719 
720  do {
721  volatile double m=0, q=0;
722  double f_sq_sum=0.0;
723 
724  for (k=0;k<lcalls_per_box;k++) {
725  double fval, bin_vol;
726 
727  random_point(x,bin,bin_vol,box,xl,xu);
728 
729  fval=func(dim,x);
730  fval*=jacbin*bin_vol;
731 
732  /* recurrence for mean and variance (sum of squares) */
733 
734  {
735  double dt=fval-m;
736  m+=dt/(k+1.0);
737  q+=dt*dt*(k/(k+1.0));
738  }
739 
740  if (mode != mode_stratified) {
741  double f_sq=fval*fval;
743  }
744  }
745 
746  intgrl+=m*lcalls_per_box;
747 
748  f_sq_sum=q*lcalls_per_box;
749 
750  tss+=f_sq_sum;
751 
752  if (mode == mode_stratified) {
753  accumulate_distribution (bin, f_sq_sum);
754  }
755 
756  } while (change_box_coord(box));
757 
758  /* Compute final results for this iteration */
759 
760  var=tss/(lcalls_per_box-1.0);
761 
762  if (var>0) {
763  wgt=1.0/var;
764  } else if (sum_wgts>0) {
765  wgt=sum_wgts/samples;
766  } else {
767  wgt=0.0;
768  }
769 
770  intgrl_sq=intgrl*intgrl;
771 
772  sig=sqrt(var);
773 
774  result=intgrl;
775  sigma=sig;
776 
777  if (wgt > 0.0) {
778  double lsum_wgts=sum_wgts;
779  double m=(sum_wgts > 0) ? (wtd_int_sum/sum_wgts) : 0;
780  double q=intgrl-m;
781 
782  samples++;
783  sum_wgts+=wgt;
784  wtd_int_sum+=intgrl*wgt;
785  chi_sum+=intgrl_sq*wgt;
786 
787  cum_int= wtd_int_sum/sum_wgts;
788  cum_sig=sqrt(1/sum_wgts);
789 
790  /* The original chisq formula from the Lepage paper is
791 
792  if ( samples > 1) {
793  chisq=(chi_sum-wtd_int_sum*cum_int)/(samples-1.0);
794  }
795 
796  This can suffer from cancellations and return a negative
797  value of chi squared. We use the new formula from
798  GSL-1.12 instead
799  */
800  if (samples==1) {
801  chisq=0;
802  } else {
803  chisq*=(samples-2.0);
804  chisq+=(wgt/(1+(wgt/lsum_wgts)))*q*q;
805  chisq/=(samples-1.0);
806  }
807 
808  } else {
809  cum_int+=(intgrl-cum_int)/(it+1.0);
810  cum_sig=0.0;
811  }
812 
813  if (this->verbose >= 1) {
814  print_res(it_num,intgrl,sig,cum_int,cum_sig,chisq);
815  if (it+1 == iterations && this->verbose > 1) {
816  print_grid(dim);
817  }
818  }
819 
820  if (this->verbose > 2) {
821  print_dist(dim);
822  }
823 
824  refine_grid ();
825 
826  if (this->verbose > 2) {
827  print_grid(dim);
828  }
829 
830  }
831 
832  /*
833  By setting stage to 1 further calls will generate independent
834  estimates based on the same grid, although it may be rebinned.
835  */
836  stage=1;
837 
838  res=cum_int;
839  err=cum_sig;
840 
841  return GSL_SUCCESS;
842  }
843 
844  virtual ~mcarlo_vegas() {}
845 
846  /// Integrate function \c func from x=a to x=b.
847  virtual int minteg_err(func_t &func, size_t ndim, const vec_t &a,
848  const vec_t &b, double &res, double &err) {
849  allocate(ndim);
850  chisq=0;
851  bins=bins_max;
852  int ret=vegas_minteg_err(0,func,ndim,a,b,res,err);
853  return ret;
854  }
855 
856  /** \brief Integrate function \c func over the hypercube from
857  \f$ x_i=a_i \f$ to \f$ x_i=b_i \f$ for
858  \f$ 0<i< \f$ ndim-1
859  */
860  virtual double minteg(func_t &func, size_t ndim, const vec_t &a,
861  const vec_t &b) {
862  double res;
863  minteg_err(func,ndim,a,b,res,this->interror);
864  return res;
865  }
866 
867  /// Return string denoting type ("mcarlo_vegas")
868  virtual const char *type() { return "mcarlo_vegas"; }
869 
870  };
871 
872 #ifndef DOXYGEN_NO_O2NS
873 }
874 #endif
875 
876 #endif
877 
virtual int vegas_minteg_err(int stage, func_t &func, size_t ndim, const vec_t &xl, const vec_t &xu, double &res, double &err)
Integrate function func from x=a to x=b.
Definition: mcarlo_vegas.h:610
ubvector_int bin
The bins for each direction.
Definition: mcarlo_vegas.h:187
virtual void print_head(unsigned long num_dim, unsigned long calls, unsigned int lit_num, unsigned int lbins, unsigned int lboxes)
Print header.
Definition: mcarlo_vegas.h:471
rng_dist_t rng_dist
The random number distribution.
Definition: mcarlo.h:68
invalid argument supplied by user
Definition: err_hnd.h:59
vec_t x
Point for function evaluation.
Definition: mcarlo_vegas.h:555
unsigned long n_points
Number of integration points (default 1000)
Definition: mcarlo.h:65
virtual double minteg(func_t &func, size_t ndim, const vec_t &a, const vec_t &b)
Integrate function func over the hypercube from to for ndim-1.
Definition: mcarlo_vegas.h:860
virtual void resize_grid(unsigned int lbins)
Resize the grid.
Definition: mcarlo_vegas.h:348
size_t dim
Number of dimensions.
Definition: mcarlo_vegas.h:170
unsigned int it_num
The total number of iterations.
Definition: mcarlo_vegas.h:207
double alpha
The stiffness of the rebinning algorithm (default 1.5)
Definition: mcarlo_vegas.h:143
ubvector_int box
The boxes for each direction.
Definition: mcarlo_vegas.h:189
unsigned int calls_per_box
Number of function calls per box.
Definition: mcarlo_vegas.h:211
virtual const char * type()
Return string denoting type ("mcarlo_vegas")
Definition: mcarlo_vegas.h:868
unsigned int iterations
Set the number of iterations (default 5)
Definition: mcarlo_vegas.h:146
virtual void init_grid(const vec_t &xl, const vec_t &xu, size_t ldim)
Initialize grid.
Definition: mcarlo_vegas.h:250
rng_t rng
The random number generator.
Definition: mcarlo.h:71
virtual void print_res(unsigned int itr, double res, double err, double cum_res, double cum_err, double chi_sq)
Print results.
Definition: mcarlo_vegas.h:488
ubvector xin
Storage for grid refinement.
Definition: mcarlo_vegas.h:178
void random_point(vec_t &lx, ubvector_int &lbin, double &bin_vol, const ubvector_int &lbox, const vec_t &xl, const vec_t &xu)
Generate a random position in a given box.
Definition: mcarlo_vegas.h:303
int verbose
Verbosity.
Definition: inte_multi.h:64
virtual void print_grid(unsigned long ldim)
Print grid.
Definition: mcarlo_vegas.h:530
#define O2SCL_ERR2(d, d2, n)
Set an error, two-string version.
Definition: err_hnd.h:281
std::string dtos(double x, int prec=6, bool auto_prec=false)
Convert a double to a string.
unsigned int bins
Number of bins.
Definition: mcarlo_vegas.h:172
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
Definition: err_hnd.h:273
Monte-Carlo integration [abstract base].
Definition: mcarlo.h:53
double vol
The volume of the current bin.
Definition: mcarlo_vegas.h:184
ubvector weight
The weight for each bin.
Definition: mcarlo_vegas.h:182
virtual int minteg_err(func_t &func, size_t ndim, const vec_t &a, const vec_t &b, double &res, double &err)
Integrate function func from x=a to x=b.
Definition: mcarlo_vegas.h:847
ubvector delx
The iteration length in each direction.
Definition: mcarlo_vegas.h:180
ubvector xi
Boundaries for each bin.
Definition: mcarlo_vegas.h:176
std::ostream * outs
The output stream to send output information (default std::cout)
Definition: mcarlo_vegas.h:160
static const size_t bins_max
Maximum number of bins.
Definition: mcarlo_vegas.h:167
double sigma
Uncertainty from last iteration.
Definition: mcarlo_vegas.h:137
virtual void reset_grid_values()
Reset grid values.
Definition: mcarlo_vegas.h:268
double interror
The uncertainty for the last integration computation.
Definition: inte_multi.h:110
virtual void refine_grid()
Refine the grid.
Definition: mcarlo_vegas.h:385
Multidimensional integration using Vegas Monte Carlo (GSL)
Definition: mcarlo_vegas.h:117
unsigned int samples
Number of samples for computing chi squared.
Definition: mcarlo_vegas.h:209
virtual int allocate(size_t ldim)
Allocate memory.
Definition: mcarlo_vegas.h:573
unsigned int boxes
Number of boxes.
Definition: mcarlo_vegas.h:174
unsigned int it_start
The starting iteration number.
Definition: mcarlo_vegas.h:205
int change_box_coord(ubvector_int &boxt)
Change box coordinates.
Definition: mcarlo_vegas.h:232
virtual void print_dist(unsigned long ldim)
Print distribution.
Definition: mcarlo_vegas.h:506
virtual void print_lim(const vec_t &xl, const vec_t &xu, unsigned long ldim)
Print limits of integration.
Definition: mcarlo_vegas.h:455
void accumulate_distribution(ubvector_int &lbin, double y)
Add the most recently generated result to the distribution.
Definition: mcarlo_vegas.h:284
double chisq
The chi-squared per degree of freedom for the weighted estimate of the integral.
Definition: mcarlo_vegas.h:157
std::function< double(size_t, const boost::numeric::ublas::vector< double > &)> multi_funct11
Multi-dimensional function typedef.
Definition: multi_funct.h:45
virtual void init_box_coord(ubvector_int &boxt)
Initialize box coordinates.
Definition: mcarlo_vegas.h:214
ubvector d
Distribution.
Definition: mcarlo_vegas.h:192
double result
Result from last iteration.
Definition: mcarlo_vegas.h:134

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).
Hosted at Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads..