gsl_vegas.h

00001 /*
00002   -------------------------------------------------------------------
00003   
00004   Copyright (C) 2006, 2007, 2008, Andrew W. Steiner
00005   
00006   This file is part of O2scl.
00007   
00008   O2scl is free software; you can redistribute it and/or modify
00009   it under the terms of the GNU General Public License as published by
00010   the Free Software Foundation; either version 3 of the License, or
00011   (at your option) any later version.
00012   
00013   O2scl is distributed in the hope that it will be useful,
00014   but WITHOUT ANY WARRANTY; without even the implied warranty of
00015   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00016   GNU General Public License for more details.
00017   
00018   You should have received a copy of the GNU General Public License
00019   along with O2scl. If not, see <http://www.gnu.org/licenses/>.
00020 
00021   -------------------------------------------------------------------
00022 */
00023 #ifndef O2SCL_GSL_VEGAS_H
00024 #define O2SCL_GSL_VEGAS_H
00025 
00026 #include <iostream>
00027 #include <o2scl/collection.h>
00028 #include <o2scl/mcarlo_inte.h>
00029 #include <o2scl/gsl_rnga.h>
00030 #include <gsl/gsl_math.h>
00031 #include <gsl/gsl_monte.h>
00032 #include <gsl/gsl_machine.h>
00033 #include <gsl/gsl_monte_vegas.h>
00034 
00035 #ifndef DOXYGENP
00036 namespace o2scl {
00037 #endif
00038   
00039   /** 
00040       \brief Multidimensional integration using Vegas Monte Carlo (GSL)
00041 
00042       The output options are a little different than the original GSL
00043       routine. The default setting of \ref mcarlo_inte::verbose is 0,
00044       which turns off all output. A verbose value of 1 prints summary
00045       information about the weighted average and final result, while a
00046       value of 2 also displays the grid coordinates. A value of 3
00047       prints information from the rebinning procedure for each
00048       iteration.
00049 
00050       Original documentation from GSL:
00051 
00052       \verbatim
00053       This is an implementation of the adaptive Monte-Carlo algorithm
00054       of G. P. Lepage, originally described in J. Comp. Phys. 27, 192(1978).
00055       The current version of the algorithm was described in the Cornell
00056       preprint CLNS-80/447 of March, 1980.
00057 
00058       This code follows most closely the c version by D.R.Yennie, coded
00059       in 1984.
00060 
00061       The input coordinates are x[j], with upper and lower limits xu[j]
00062       and xl[j].  The integration length in the j-th direction is
00063       delx[j].  Each coordinate x[j] is rescaled to a variable y[j] in
00064       the range 0 to 1.  The range is divided into bins with boundaries
00065       xi[i][j], where i=0 corresponds to y=0 and i=bins to y=1.  The grid
00066       is refined (ie, bins are adjusted) using d[i][j] which is some
00067       variation on the squared sum.  A third parameter used in defining
00068       the real coordinate using random numbers is called z.  It ranges
00069       from 0 to bins.  Its integer part gives the lower index of the bin
00070       into which a call is to be placed, and the remainder gives the
00071       location inside the bin.
00072 
00073       When stratified sampling is used the bins are grouped into boxes,
00074       and the algorithm allocates an equal number of function calls to
00075       each box.
00076 
00077       The variable alpha controls how "stiff" the rebinning algorithm is.  
00078       alpha = 0 means never change the grid.  Alpha is typically set between
00079       1 and 2.
00080       \endverbatim
00081 
00082       \todo Need to double check that the verbose output is good
00083       for all settings of verbose.
00084       \todo BINS_MAX and bins_max are somehow duplicates. Fix this.
00085 
00086   */
00087   template<class param_t, class func_t=multi_funct<param_t>, 
00088     class rng_t=gsl_rnga, class vec_t=ovector_view, class alloc_vec_t=ovector,
00089     class alloc_t=ovector_alloc> class gsl_vegas : 
00090   public mcarlo_inte<param_t,func_t,rng_t,vec_t> {
00091 
00092     public:
00093     
00094     /// \name Integration mode (default is mode_importance)
00095     //@{
00096     int mode;
00097     static const int mode_importance = 1;
00098     static const int mode_importance_only = 0;
00099     static const int mode_stratified = -1;
00100     //@}
00101 
00102     /// Result from last iteration
00103     double result;
00104 
00105     /// Uncertainty from last iteration
00106     double sigma;
00107 
00108     /** 
00109         \brief The stiffness of the rebinning algorithm (default 1.5)
00110 
00111         This usual range is between 1 and 2.
00112     */
00113     double alpha;
00114 
00115     /// Set the number of iterations (default 5)
00116     unsigned int iterations;
00117 
00118     /**
00119        \brief The chi-squared per degree of freedom for the weighted
00120        estimate of the integral
00121 
00122        From GSL documentation:
00123        \verbatim
00124        The value of CHISQ should be close to 1.  A value of CHISQ
00125        which differs significantly from 1 indicates that the values
00126        from different iterations are inconsistent.  In this case the
00127        weighted error will be under-estimated, and further iterations
00128        of the algorithm are needed to obtain reliable results.
00129        \endverbatim
00130     */
00131     double chisq;
00132     
00133     /// The output stream to send output information (default \c std::cout)
00134     std::ostream *outs;
00135 
00136 #ifndef DOXYGEN_INTERNAL
00137 
00138     protected:
00139 
00140     /// Desc
00141     typedef int coord;
00142 
00143     /// Desc
00144     static const int BINS_MAX=50;
00145     
00146     /* grid */
00147     size_t dim;
00148     size_t bins_max;
00149     unsigned int bins;
00150     /* these are both counted along the axes */
00151     unsigned int boxes; 
00152     double *xi;
00153     double *xin;
00154     double *delx;
00155     double *weight;
00156     double vol;
00157 
00158     int *bin;
00159     int *box;
00160 
00161     /* distribution */
00162     double *d;
00163 
00164     /* scratch variables preserved between calls to vegas1/2/3  */
00165     double jac;
00166     double wtd_int_sum;
00167     double sum_wgts;
00168     double chi_sum;
00169     
00170     unsigned int it_start;
00171     unsigned int it_num;
00172     unsigned int samples;
00173     unsigned int calls_per_box;
00174     
00175     /// Initialize box coordinates
00176     virtual void init_box_coord(coord boxt[]) {
00177       size_t i;
00178       for (i=0;i<dim;i++) {
00179         boxt[i]=0;
00180       }
00181       return;
00182     }
00183   
00184     /**
00185        \brief Change box coordinates
00186 
00187        Steps through the box coord like
00188        {0,0}, {0, 1}, {0, 2}, {0, 3}, {1, 0}, {1, 1}, {1, 2}, ...
00189 
00190        This is among the member functions that is not virtual
00191        because it is part of the innermost loop.
00192     */
00193     int change_box_coord(coord boxt[]) {
00194       int j = dim - 1;
00195       int ng = boxes;
00196       
00197       while (j >= 0) {
00198 
00199         boxt[j] = (boxt[j] + 1) % ng;
00200         if (boxt[j] != 0) {
00201           return 1;
00202         }
00203         j--;
00204 
00205       }
00206       
00207       return 0;
00208     }
00209  
00210     /// Initialize grid
00211     virtual void init_grid(const vec_t &xl, const vec_t &xu, size_t ldim) {
00212       size_t j;
00213       vol = 1.0;
00214       bins = 1;
00215       
00216       for (j = 0; j < ldim; j++) {
00217         double dx = xu[j] - xl[j];
00218         delx[j] = dx;
00219         vol *= dx;
00220         
00221         xi[j]=0.0;
00222         xi[dim+j]=1.0;
00223       }
00224       
00225       return;
00226     }
00227 
00228     /// Reset grid values
00229     virtual void reset_grid_values() {
00230       size_t i, j;
00231       
00232       for (i = 0; i < bins; i++) {
00233         for (j = 0; j < dim; j++) {
00234           d[i*dim+j] = 0.0;
00235         }
00236       }
00237       return;
00238     }
00239 
00240     /** 
00241         \brief Add the most recently generated result to the distribution
00242 
00243         This is among the member functions that is not virtual
00244         because it is part of the innermost loop.
00245     */
00246     void accumulate_distribution(coord lbin[], double y) {
00247       size_t j;
00248       
00249       for (j=0;j<dim;j++) {
00250         int i=lbin[j];
00251         d[i*dim+j]+=y;
00252       }
00253       return;
00254     }
00255     
00256     /**
00257        \brief Generate a random position in a given box
00258        
00259        Use the random number generator r to return a random position x
00260        in a given box.  The value of bin gives the bin location of the
00261        random position (there may be several bins within a given box)
00262        
00263        This is among the member functions that is not virtual
00264        because it is part of the innermost loop.
00265     */
00266     void random_point(vec_t &lx, coord lbin[], double *bin_vol,
00267                       const coord lbox[], const vec_t &xl, 
00268                       const vec_t &xu) {
00269 
00270       double lvol = 1.0;
00271 
00272       size_t j;
00273 
00274       for (j = 0; j < dim; ++j) {
00275 
00276         // The equivalent of gsl_rng_uniform_pos()
00277         double rdn;
00278         do { rdn=this->def_rng.random(); } while (rdn==0);
00279         
00280         /* lbox[j] + ran gives the position in the box units, while z
00281            is the position in bin units.  */
00282         double z = ((lbox[j] + rdn) / boxes) * bins;
00283         
00284         int k = (int)z;
00285         
00286         double y, bin_width;
00287         
00288         lbin[j] = k;
00289 
00290         if (k == 0) {
00291           bin_width = xi[dim+j];
00292           y = z * bin_width;
00293         } else {
00294           bin_width = xi[(k+1)*dim+j]-xi[k*dim+j];
00295           y = xi[k*dim+j] + (z - k) * bin_width;
00296         }
00297         
00298         lx[j] = xl[j] + y * delx[j];
00299         
00300         lvol *= bin_width;
00301       }
00302       
00303       *bin_vol = lvol;
00304 
00305       return;
00306     }
00307 
00308     /// Resize the grid
00309     virtual void resize_grid(unsigned int lbins) {
00310       size_t j, k;
00311 
00312       /* weight is ratio of bin sizes */
00313 
00314       double pts_per_bin = (double) bins / (double) lbins;
00315 
00316       for (j = 0; j < dim; j++) {
00317         double xold;
00318         double xnew = 0;
00319         double dw = 0;
00320         int i = 1;
00321 
00322         for (k = 1; k <= bins; k++) {
00323           dw += 1.0;
00324           xold = xnew;
00325           xnew = xi[k*dim+j];
00326 
00327           for (; dw > pts_per_bin; i++) {
00328             dw -= pts_per_bin;
00329             (xin[i]) = xnew - (xnew - xold) * dw;
00330           }
00331         }
00332 
00333         for (k = 1 ; k < lbins; k++) {
00334           xi[k*dim+j] = xin[k];
00335         }
00336 
00337         xi[lbins*dim+j] = 1;
00338       }
00339       
00340       bins=lbins;
00341 
00342       return;
00343     }
00344 
00345     /// Refine the grid
00346     virtual void refine_grid() {
00347       size_t i, j, k;
00348 
00349       for (j = 0; j < dim; j++) {
00350         double grid_tot_j, tot_weight;
00351 
00352         double oldg = d[j];
00353         double newg = d[dim+j];
00354 
00355         d[j] = (oldg + newg) / 2;
00356         grid_tot_j = d[j];
00357 
00358         /* This implements gs[i][j] = (gs[i-1][j]+gs[i][j]+gs[i+1][j])/3 */
00359 
00360         for (i = 1; i < bins - 1; i++) {
00361           double rc = oldg + newg;
00362           oldg = newg;
00363           newg = d[(i+1)*dim+j];
00364           d[i*dim+j] = (rc + newg) / 3;
00365           grid_tot_j += d[i*dim+j];
00366         }
00367         d[(bins-1)*dim+j] = (newg + oldg) / 2;
00368 
00369         grid_tot_j += d[(bins-1)*dim+j];
00370 
00371         tot_weight = 0;
00372 
00373         for (i = 0; i < bins; i++) {
00374           weight[i] = 0;
00375 
00376           if (d[i*dim+j] > 0) {
00377             oldg = grid_tot_j / d[i*dim+j];
00378             /* damped change */
00379             weight[i] = pow (((oldg - 1) / oldg / log (oldg)), alpha);
00380           }
00381               
00382           tot_weight += weight[i];
00383         }
00384 
00385         {
00386           double pts_per_bin = tot_weight / bins;
00387 
00388           double xold;
00389           double xnew = 0;
00390           double dw = 0;
00391           i = 1;
00392 
00393           for (k = 0; k < bins; k++) {
00394             dw += weight[k];
00395             xold = xnew;
00396             xnew = xi[(k+1)*dim+j];
00397 
00398             for (; dw > pts_per_bin; i++) {
00399               dw -= pts_per_bin;
00400               (xin[i]) = xnew - (xnew - xold) * dw / weight[k];
00401             }
00402           }
00403 
00404           for (k = 1 ; k < bins ; k++) {
00405             xi[k*dim+j] = xin[k];
00406           }
00407             
00408           xi[bins*dim+j]=1;
00409         }
00410       }
00411       return;
00412     }
00413 
00414     /// Print limits of integration
00415     virtual void print_lim(const vec_t &xl, const vec_t &xu, 
00416                            unsigned long ldim) {
00417                    
00418       unsigned long j;
00419       
00420       (*outs) << "The limits of integration are:\n" << std::endl;
00421       for (j = 0; j < ldim; ++j) {
00422         (*outs) << "xl[" << j << "]=" << xl[j] << "    xu[" << j 
00423                 <<  "]=" << xu[j] << std::endl;
00424       }
00425       (*outs) << std::endl;
00426 
00427       return;
00428     }
00429 
00430     /// Print header
00431     virtual void print_head(unsigned long num_dim, unsigned long calls,
00432                             unsigned int lit_num, unsigned int lbins, 
00433                             unsigned int lboxes)  {
00434       (*outs) << "num_dim=" << num_dim << ", calls=" << calls 
00435       << ", it_num=" << lit_num << ", max_it_num=" 
00436       << iterations << ", verb=" << this->verbose << ", alph=" << alpha 
00437       << ",\n mode=" << mode << ", boxes=" << lboxes
00438       << "\n\n    single.......iteration               "
00439       << "accumulated......results\n"
00440       << "iter   integral     sigma            integral   "
00441       << "  sigma        chi-sq/it" << std::endl;
00442       
00443       return;
00444     }
00445     
00446     /// Print results
00447     virtual void print_res(unsigned int itr, double res, 
00448                            double err, double cum_res, double cum_err,
00449                            double chi_sq) {
00450       (*outs).precision(5);
00451       (*outs) << itr << "     ";
00452       outs->setf(std::ios::showpos);
00453       (*outs) << res << "  ";
00454       outs->unsetf(std::ios::showpos);
00455       (*outs) << err << "     "; 
00456       outs->setf(std::ios::showpos);
00457       (*outs) << cum_res << "  ";
00458       outs->unsetf(std::ios::showpos);
00459       (*outs) << cum_err << "  " << chi_sq << std::endl;
00460       (*outs).precision(6);
00461       return;
00462     }
00463  
00464     /// Print distribution
00465     virtual void print_dist(unsigned long ldim) {
00466       unsigned long i, j;
00467       int p = this->verbose;
00468       if (p < 2) return;
00469       
00470       for (j = 0; j < ldim; ++j) {
00471         (*outs) << "\n axis " << j << std::endl;
00472         (*outs) << "       x     g" << std::endl;
00473         for (i = 0; i < bins; i++) {
00474           (*outs) << "weight [" << (xi[(i)*dim + (j)]) << " , " 
00475                   << xi[(i+1)*dim+j] << "] = ";
00476           (*outs) << " " << (d[(i)*dim + (j)]) << std::endl;
00477           
00478         }
00479         (*outs) << std::endl;
00480       }
00481       (*outs) << std::endl;
00482       return;
00483     }
00484  
00485     /// Print grid
00486     virtual void print_grid(unsigned long ldim) {
00487       unsigned long i, j;
00488       int p = this->verbose;
00489       if (p < 2) {
00490         return;
00491       }
00492       
00493       for (j = 0; j < ldim; ++j) {
00494         (*outs) << "\n axis " << j << std::endl;
00495         (*outs) << "      x     " << std::endl;
00496         for (i = 0; i <= bins; i++)      {
00497           (*outs) << (xi[(i)*dim + (j)]) << " ";
00498           if (i % 5 == 4) {
00499             (*outs) << std::endl;
00500           }
00501         }
00502         (*outs) << std::endl;
00503       }
00504       (*outs) << std::endl;
00505       return;
00506     }
00507 
00508     /// Memory allocation object
00509     alloc_t ao;
00510 
00511     /// Point for function evaluation
00512     alloc_vec_t x;
00513 
00514 #endif
00515 
00516     public:
00517 
00518     gsl_vegas() {
00519       this->verbose = 0;
00520       outs=&std::cout;
00521       
00522       alpha=1.5;
00523       iterations=5;
00524       mode=mode_importance;
00525       chisq=0;
00526       bins=bins_max;
00527     }
00528 
00529     /// Allocate memory
00530     virtual int allocate(size_t ldim) {
00531     
00532       delx=(double *)malloc(ldim*sizeof(double));
00533       if (delx == 0) {
00534         GSL_ERROR_VAL("failed to allocate space for delx", GSL_ENOMEM, 0);
00535       }
00536 
00537       d=(double *)malloc(BINS_MAX*ldim*sizeof(double));
00538       if (d == 0) {
00539         std::free(delx);
00540         GSL_ERROR_VAL("failed to allocate space for d", GSL_ENOMEM, 0);
00541       }
00542 
00543       xi=(double *)malloc((BINS_MAX+1)*ldim*sizeof(double));
00544       if (xi == 0) {
00545         std::free(d);
00546         std::free(delx);
00547         GSL_ERROR_VAL("failed to allocate space for xi", GSL_ENOMEM, 0);
00548       }
00549 
00550       xin=(double *)malloc((BINS_MAX+1)*sizeof (double));
00551       if (xin == 0) {
00552         std::free(xi);
00553         std::free(d);
00554         std::free(delx);
00555         GSL_ERROR_VAL("failed to allocate space for xin", GSL_ENOMEM, 0);
00556       }
00557 
00558       weight=(double *)malloc(BINS_MAX*sizeof (double));
00559       if (weight == 0) {
00560         std::free(xin);
00561         std::free(xi);
00562         std::free(d);
00563         std::free(delx);
00564         GSL_ERROR_VAL("failed to allocate space for xin", GSL_ENOMEM, 0);
00565       }
00566 
00567       box=(coord *)malloc(ldim*sizeof (coord));
00568       if (box == 0) {
00569         std::free(weight);
00570         std::free(xin);
00571         std::free(xi);
00572         std::free(d);
00573         std::free(delx);
00574         GSL_ERROR_VAL("failed to allocate space for box", GSL_ENOMEM, 0);
00575       }
00576 
00577       bin=(coord *)malloc(ldim*sizeof (coord));
00578       if (bin == 0) {
00579         std::free(box);
00580         std::free(weight);
00581         std::free(xin);
00582         std::free(xi);
00583         std::free(d);
00584         std::free(delx);
00585         GSL_ERROR_VAL("failed to allocate space for bin", GSL_ENOMEM, 0);
00586       }
00587 
00588       ao.allocate(x,ldim);
00589 
00590       dim=ldim;
00591       bins_max=BINS_MAX;
00592 
00593       return 0;
00594     }
00595 
00596     /// Free allocated memory
00597     virtual int free() {
00598       std::free(delx);
00599       std::free(d);
00600       std::free(xi);
00601       std::free(xin);
00602       std::free(weight);
00603       std::free(box);
00604       std::free(bin);
00605       ao.free(x);
00606       return 0;
00607     }
00608 
00609     /** 
00610         \brief Integrate function \c func from x=a to x=b.
00611 
00612         Original documentation from GSL:
00613         
00614         Normally, `stage = 0' which begins with a new uniform grid and
00615         empty weighted average.  Calling vegas with `stage = 1'
00616         retains the grid from the previous run but discards the
00617         weighted average, so that one can "tune" the grid using a
00618         relatively small number of points and then do a large run with
00619         `stage = 1' on the optimized grid.  Setting `stage = 2' keeps
00620         the grid and the weighted average from the previous run, but
00621         may increase (or decrease) the number of histogram bins in the
00622         grid depending on the number of calls available.  Choosing
00623         `stage = 3' enters at the main loop, so that nothing is
00624         changed, and is equivalent to performing additional iterations
00625         in a previous call.
00626     */
00627     virtual int vegas_minteg_err(int stage, func_t &func, size_t ndim, 
00628                                  const vec_t &xl, const vec_t &xu, param_t &pa,
00629                                  double &res, double &err) {
00630 
00631       size_t calls=this->n_points;
00632 
00633       double cum_int, cum_sig;
00634       size_t i, k, it;
00635         
00636       for (i = 0; i < dim; i++) {
00637         if (xu[i] <= xl[i]) {
00638           GSL_ERROR ("xu must be greater than xl", GSL_EINVAL);
00639         }
00640 
00641         if (xu[i] - xl[i] > GSL_DBL_MAX) {
00642           GSL_ERROR ("Range of integration is too large, please rescale",
00643                      GSL_EINVAL);
00644         }
00645       }
00646 
00647       if (stage == 0) {
00648         init_grid(xl,xu,dim);
00649         if (this->verbose>=1) {
00650           print_lim(xl,xu,dim);
00651         }
00652       }
00653       
00654       if (stage<=1) {
00655         wtd_int_sum = 0;
00656         sum_wgts = 0;
00657         chi_sum = 0;
00658         it_num = 1;
00659         samples = 0;
00660       }
00661       
00662       if (stage <= 2) {
00663 
00664         unsigned int lbins = bins_max;
00665         unsigned int lboxes = 1;
00666 
00667         if (mode != mode_importance_only) {
00668 
00669           /* shooting for 2 calls/box */
00670           
00671           // The original GSL code was:
00672           // boxes = floor (pow (calls / 2.0, 1.0 / dim));
00673           // but floor returns double on my machine, so 
00674           // we explicitly typecast here
00675           
00676           lboxes = ((unsigned int)(floor(pow(calls/2.0,1.0/dim))));
00677           mode = mode_importance;
00678 
00679           if (2 * lboxes >=  bins_max) {
00680             /* if bins/box < 2 */
00681             int box_per_bin = GSL_MAX(lboxes/bins_max,1);
00682 
00683             lbins = GSL_MIN(lboxes/box_per_bin,bins_max);
00684             lboxes = box_per_bin*lbins;
00685 
00686             mode = mode_stratified;
00687           }
00688 
00689         }
00690 
00691         double tot_boxes = pow ((double) lboxes, (double) dim);
00692         calls_per_box=((unsigned int)(GSL_MAX(calls/tot_boxes,2)));
00693         calls = ((size_t)( calls_per_box * tot_boxes));
00694         
00695         /* total volume of x-space/(avg num of calls/bin) */
00696         jac =  vol * pow ((double) lbins, (double) dim) / calls;
00697          
00698         boxes=lboxes;
00699 
00700         /* If the number of bins changes from the previous invocation, bins
00701            are expanded or contracted accordingly, while preserving bin
00702            density */
00703         
00704         if (lbins!=bins) {
00705           resize_grid(lbins);
00706           if (this->verbose > 2) print_grid(dim);
00707         }
00708         if (this->verbose >= 1) {
00709           print_head(dim,calls,it_num,bins,boxes);
00710         }
00711       }
00712 
00713       it_start = it_num;
00714 
00715       cum_int = 0.0;
00716       cum_sig = 0.0;
00717 
00718       chisq = 0.0;
00719       
00720       for (it = 0; it < iterations; it++) {
00721 
00722         double intgrl = 0.0, intgrl_sq = 0.0;
00723         double sig = 0.0;
00724         double wgt;
00725         size_t lcalls_per_box =  calls_per_box;
00726         double jacbin =  jac;
00727 
00728         it_num = it_start + it;
00729 
00730         reset_grid_values();
00731         init_box_coord(box);
00732 
00733         do {
00734           double m = 0, q = 0;
00735           double f_sq_sum = 0.0;
00736 
00737           for (k = 0; k < lcalls_per_box; k++) {
00738             double fval, bin_vol;
00739             
00740             random_point(x,bin,&bin_vol,box,xl,xu);
00741             
00742             func(dim,x,fval,pa);
00743             fval*=jacbin*bin_vol;
00744 
00745             /* recurrence for mean and variance */
00746 
00747             {
00748               double dt = fval - m;
00749               m += dt / (k + 1.0);
00750               q += dt * dt * (k / (k + 1.0));
00751             }
00752 
00753             if ( mode != mode_stratified) {
00754               double f_sq = fval * fval;
00755               accumulate_distribution (bin, f_sq);
00756             }
00757           }
00758 
00759           intgrl += m * lcalls_per_box;
00760 
00761           f_sq_sum = q * lcalls_per_box ;
00762 
00763           sig += f_sq_sum ;
00764 
00765           if ( mode == mode_stratified) {
00766             accumulate_distribution (bin, f_sq_sum);
00767           }
00768 
00769         } while (change_box_coord(box));
00770 
00771         /* Compute final results for this iteration   */
00772         
00773         sig=sig/(lcalls_per_box-1.0);
00774         
00775         if (sig > 0) {
00776           wgt = 1.0 / sig;
00777         } else if ( sum_wgts > 0) {
00778           wgt =  sum_wgts /  samples;
00779         } else {
00780           wgt = 0.0;
00781         }
00782         
00783         intgrl_sq = intgrl * intgrl;
00784         
00785         result = intgrl;
00786         sigma = sqrt(sig);
00787         
00788         if (wgt > 0.0) {
00789           samples++ ;
00790           sum_wgts += wgt;
00791           wtd_int_sum += intgrl * wgt;
00792           chi_sum += intgrl_sq * wgt;
00793 
00794           cum_int =  wtd_int_sum /  sum_wgts;
00795           cum_sig = sqrt (1 /  sum_wgts);
00796 
00797           if ( samples > 1) {
00798             chisq=(chi_sum-wtd_int_sum*cum_int)/(samples-1.0);
00799           }
00800         } else {
00801           cum_int += (intgrl - cum_int) / (it + 1.0);
00802           cum_sig = 0.0;
00803         }         
00804         
00805         if ( this->verbose >= 1) {
00806           print_res(it_num, intgrl, sqrt (sig), 
00807                     cum_int, cum_sig,  chisq);
00808           if (it + 1 ==  iterations &&  this->verbose > 1) {
00809             print_grid(dim);
00810           }
00811         }
00812 
00813         if ( this->verbose > 2) {
00814           print_dist(dim);
00815         }
00816 
00817         refine_grid ();
00818 
00819         if ( this->verbose > 2) {
00820           print_grid(dim);
00821         }
00822 
00823       }
00824 
00825       /* 
00826          By setting stage to 1 further calls will generate independent
00827          estimates based on the same grid, although it may be rebinned. 
00828       */
00829 
00830       stage = 1;  
00831 
00832       res = cum_int;
00833       err = cum_sig;
00834 
00835       return GSL_SUCCESS;
00836     }
00837 
00838     virtual ~gsl_vegas() {}
00839   
00840     /// Integrate function \c func from x=a to x=b.
00841     virtual int minteg_err(func_t &func, size_t ndim, const vec_t &a, 
00842                            const vec_t &b, param_t &pa,
00843                            double &res, double &err) {
00844       allocate(ndim);
00845       chisq = 0;
00846       bins = bins_max;
00847       int ret=vegas_minteg_err(0,func,ndim,a,b,pa,res,err);
00848       free();
00849       return ret;
00850     }
00851       
00852     /// Return string denoting type ("gsl_vegas")
00853     virtual const char *type() { return "gsl_vegas"; }
00854       
00855   };
00856 
00857 #ifndef DOXYGENP
00858 }
00859 #endif
00860 
00861 #endif
00862 

Documentation generated with Doxygen and provided under the GNU Free Documentation License. See License Information for details.

Project hosting provided by SourceForge.net Logo, O2scl Sourceforge Project Page