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