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