00001 /* 00002 ------------------------------------------------------------------- 00003 00004 Copyright (C) 2006, 2007, 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_VEGAS2_H 00024 #define O2SCL_GSL_VEGAS2_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 plain Vegas Carlo (GSL) 00041 00042 (Upgrade of gsl_vegas. Not working yet.) 00043 00044 Original documentation from GSL: 00045 00046 This is an implementation of the adaptive Monte-Carlo algorithm 00047 of G. P. Lepage, originally described in J. Comp. Phys. 27, 192(1978). 00048 The current version of the algorithm was described in the Cornell 00049 preprint CLNS-80/447 of March, 1980. 00050 00051 This code follows most closely the c version by D.R.Yennie, coded 00052 in 1984. 00053 00054 The input coordinates are x[j], with upper and lower limits xu[j] 00055 and xl[j]. The integration length in the j-th direction is 00056 delx[j]. Each coordinate x[j] is rescaled to a variable y[j] in 00057 the range 0 to 1. The range is divided into bins with boundaries 00058 xi[i][j], where i=0 corresponds to y=0 and i=bins to y=1. The grid 00059 is refined (ie, bins are adjusted) using d[i][j] which is some 00060 variation on the squared sum. A third parameter used in defining 00061 the real coordinate using random numbers is called z. It ranges 00062 from 0 to bins. Its integer part gives the lower index of the bin 00063 into which a call is to be placed, and the remainder gives the 00064 location inside the bin. 00065 00066 When stratified sampling is used the bins are grouped into boxes, 00067 and the algorithm allocates an equal number of function calls to 00068 each box. 00069 00070 The variable alpha controls how "stiff" the rebinning algorithm is. 00071 alpha = 0 means never change the grid. Alpha is typically set between 00072 1 and 2. 00073 00074 */ 00075 template<class param_t, class func_t=multi_funct<param_t>, 00076 class rng_t=gsl_rnga, class vec_t=ovector_view, class alloc_vec_t=ovector, 00077 class alloc_t=ovector_alloc> class gsl_vegas2 : 00078 public mcarlo_inte<param_t,func_t,rng_t,vec_t> { 00079 00080 public: 00081 00082 static const int GSL_VEGAS2_MODE_IMPORTANCE = 1; 00083 static const int GSL_VEGAS2_MODE_IMPORTANCE_ONLY = 0; 00084 static const int GSL_VEGAS2_MODE_STRATIFIED = -1; 00085 00086 protected: 00087 00088 typedef int coord; 00089 00090 static const int BINS_MAX=50; 00091 00092 /* grid */ 00093 size_t dim; 00094 size_t bins_max; 00095 unsigned int bins; 00096 unsigned int boxes; /* these are both counted along the axes */ 00097 double * xi; 00098 double * xin; 00099 double * delx; 00100 double * weight; 00101 double vol; 00102 00103 int * bin; 00104 int * box; 00105 00106 /* distribution */ 00107 double * d; 00108 00109 /* control variables */ 00110 double alpha; 00111 int mode; 00112 unsigned int iterations; 00113 int stage; 00114 00115 /* scratch variables preserved between calls to vegas1/2/3 */ 00116 double jac; 00117 double wtd_int_sum; 00118 double sum_wgts; 00119 double chi_sum; 00120 double chisq; 00121 00122 double result; 00123 double sigma; 00124 00125 unsigned int it_start; 00126 unsigned int it_num; 00127 unsigned int samples; 00128 unsigned int calls_per_box; 00129 00130 public: 00131 00132 /// Desc 00133 int allocate(size_t ldim) { 00134 00135 delx = (double *) malloc (ldim * sizeof (double)); 00136 00137 if (delx == 0) { 00138 GSL_ERROR_VAL ("failed to allocate space for delx", GSL_ENOMEM, 0); 00139 } 00140 00141 d = (double *) malloc (BINS_MAX * ldim * sizeof (double)); 00142 00143 if (d == 0) { 00144 std::free(delx); 00145 GSL_ERROR_VAL ("failed to allocate space for d", GSL_ENOMEM, 0); 00146 } 00147 00148 xi = (double *) malloc ((BINS_MAX + 1) * ldim * sizeof (double)); 00149 00150 if (xi == 0) { 00151 std::free(d); 00152 std::free(delx); 00153 GSL_ERROR_VAL ("failed to allocate space for xi", GSL_ENOMEM, 0); 00154 } 00155 00156 xin = (double *) malloc ((BINS_MAX + 1) * sizeof (double)); 00157 00158 if (xin == 0) { 00159 std::free(xi); 00160 std::free(d); 00161 std::free(delx); 00162 GSL_ERROR_VAL ("failed to allocate space for xin", GSL_ENOMEM, 0); 00163 } 00164 00165 weight = (double *) malloc (BINS_MAX * sizeof (double)); 00166 00167 if (weight == 0) { 00168 std::free(xin); 00169 std::free(xi); 00170 std::free(d); 00171 std::free(delx); 00172 GSL_ERROR_VAL ("failed to allocate space for xin", GSL_ENOMEM, 0); 00173 } 00174 00175 box = (coord *) malloc (ldim * sizeof (coord)); 00176 00177 if (box == 0) { 00178 std::free(weight); 00179 std::free(xin); 00180 std::free(xi); 00181 std::free(d); 00182 std::free(delx); 00183 GSL_ERROR_VAL ("failed to allocate space for box", GSL_ENOMEM, 0); 00184 } 00185 00186 bin = (coord *) malloc (ldim * sizeof (coord)); 00187 00188 if (bin == 0) { 00189 std::free(box); 00190 std::free(weight); 00191 std::free(xin); 00192 std::free(xi); 00193 std::free(d); 00194 std::free(delx); 00195 GSL_ERROR_VAL ("failed to allocate space for bin", GSL_ENOMEM, 0); 00196 } 00197 00198 dim = ldim; 00199 bins_max = BINS_MAX; 00200 00201 init(); 00202 00203 return 0; 00204 } 00205 00206 /* Set some default values */ 00207 int init() { 00208 stage = 0; 00209 alpha = 1.5; 00210 this->verbose = -1; 00211 iterations = 5; 00212 mode = GSL_VEGAS2_MODE_IMPORTANCE; 00213 chisq = 0; 00214 bins = bins_max; 00215 00216 return GSL_SUCCESS; 00217 } 00218 00219 /// Desc 00220 void free() { 00221 std::free(delx); 00222 std::free(d); 00223 std::free(xi); 00224 std::free(xin); 00225 std::free(weight); 00226 std::free(box); 00227 std::free(bin); 00228 } 00229 00230 protected: 00231 00232 /// Desc 00233 void init_box_coord(coord boxt[]) { 00234 size_t i; 00235 00236 for (i = 0; i < dim; i++){ 00237 boxt[i] = 0; 00238 } 00239 return; 00240 } 00241 00242 /** 00243 \brief Desc 00244 00245 change_box_coord steps through the box coord like 00246 {0,0}, {0, 1}, {0, 2}, {0, 3}, {1, 0}, {1, 1}, {1, 2}, ... 00247 */ 00248 int change_box_coord(coord boxt[]) { 00249 int j = dim - 1; 00250 int ng = boxes; 00251 00252 while (j >= 0) { 00253 boxt[j] = (boxt[j] + 1) % ng; 00254 00255 if (boxt[j] != 0) { 00256 return 1; 00257 } 00258 00259 j--; 00260 } 00261 00262 return 0; 00263 } 00264 00265 /// Desc 00266 void init_grid(const vec_t &xl, const vec_t &xu, size_t ldim) { 00267 size_t j; 00268 vol = 1.0; 00269 bins = 1; 00270 00271 for (j = 0; j < ldim; j++) { 00272 double dx = xu[j] - xl[j]; 00273 delx[j] = dx; 00274 vol *= dx; 00275 00276 (xi[(j)])=0.0; 00277 (xi[(1)*dim+(j)])=0.0; 00278 } 00279 00280 std::cout << "vol: " << delx[0] << " " << vol << std::endl; 00281 00282 return; 00283 } 00284 00285 /// Desc 00286 void reset_grid_values() { 00287 size_t i, j; 00288 00289 for (i = 0; i < bins; i++) { 00290 for (j = 0; j < dim; j++) { 00291 d[i*dim+j] = 0.0; 00292 } 00293 } 00294 return; 00295 } 00296 00297 /// Desc 00298 void accumulate_distribution(coord lbin[], double y) { 00299 size_t j; 00300 00301 for (j = 0; j < dim; j++){ 00302 int i = lbin[j]; 00303 d[i*dim+j] += y; 00304 } 00305 return; 00306 } 00307 00308 /// Desc 00309 void random_point(vec_t &lx, coord lbin[], double *bin_vol, 00310 const coord lbox[], const vec_t &xl, 00311 const vec_t &xu, rng_t &r) { 00312 00313 /* Use the random number generator r to return a random position x 00314 in a given box. The value of bin gives the bin location of the 00315 random position (there may be several bins within a given box) */ 00316 00317 double lvol = 1.0; 00318 00319 size_t j; 00320 00321 for (j = 0; j < dim; ++j) { 00322 /* lbox[j] + ran gives the position in the box units, while z 00323 is the position in bin units. */ 00324 00325 double rdn; 00326 do { rdn=r.random(); } while (rdn==0); 00327 double z = ((lbox[j] + rdn) / boxes) * bins; 00328 00329 int k = (int)z; 00330 00331 double y, bin_width; 00332 00333 lbin[j] = k; 00334 00335 if (k == 0) { 00336 bin_width = xi[dim+j]; 00337 y = z * bin_width; 00338 } 00339 else 00340 { 00341 bin_width = xi[(k+1)*dim+j]-xi[k*dim+j]; 00342 y = xi[k*dim+j] + (z - k) * bin_width; 00343 } 00344 00345 lx[j] = xl[j] + y * delx[j]; 00346 00347 lvol *= bin_width; 00348 } 00349 00350 *bin_vol = lvol; 00351 return; 00352 } 00353 00354 /// Desc 00355 void resize_grid(unsigned int lbins) { 00356 size_t j, k; 00357 00358 /* weight is ratio of bin sizes */ 00359 00360 double pts_per_bin = (double) bins / (double) lbins; 00361 00362 for (j = 0; j < dim; j++) { 00363 double xold; 00364 double xnew = 0; 00365 double dw = 0; 00366 int i = 1; 00367 00368 for (k = 1; k <= bins; k++) { 00369 dw += 1.0; 00370 xold = xnew; 00371 xnew = xi[k*dim+j]; 00372 00373 for (; dw > pts_per_bin; i++) { 00374 dw -= pts_per_bin; 00375 (xin[i]) = xnew - (xnew - xold) * dw; 00376 } 00377 } 00378 00379 for (k = 1 ; k < lbins; k++) { 00380 xi[k*dim+j] = xin[k]; 00381 } 00382 00383 xi[lbins*dim+j] = 1; 00384 } 00385 00386 return; 00387 } 00388 00389 /// Desc 00390 void refine_grid() { 00391 size_t i, j, k; 00392 00393 for (j = 0; j < dim; j++) { 00394 double grid_tot_j, tot_weight; 00395 00396 double oldg = d[j]; 00397 double newg = d[dim+j]; 00398 00399 d[j] = (oldg + newg) / 2; 00400 grid_tot_j = d[j]; 00401 00402 /* This implements gs[i][j] = (gs[i-1][j]+gs[i][j]+gs[i+1][j])/3 */ 00403 00404 for (i = 1; i < bins - 1; i++) { 00405 double rc = oldg + newg; 00406 oldg = newg; 00407 newg = d[(i+1)*dim+j]; 00408 d[i*dim+j] = (rc + newg) / 3; 00409 grid_tot_j += d[i*dim+j]; 00410 } 00411 d[(bins-1)*dim+j] = (newg + oldg) / 2; 00412 00413 grid_tot_j += d[(bins-1)*dim+j]; 00414 00415 tot_weight = 0; 00416 00417 for (i = 0; i < bins; i++) { 00418 weight[i] = 0; 00419 00420 if (d[i*dim+j] > 0) { 00421 oldg = grid_tot_j / d[i*dim+j]; 00422 /* damped change */ 00423 weight[i] = pow (((oldg - 1) / oldg / log (oldg)), alpha); 00424 } 00425 00426 tot_weight += weight[i]; 00427 } 00428 00429 { 00430 double pts_per_bin = tot_weight / bins; 00431 00432 double xold; 00433 double xnew = 0; 00434 double dw = 0; 00435 i = 1; 00436 00437 for (k = 0; k < bins; k++) { 00438 dw += weight[k]; 00439 xold = xnew; 00440 xnew = xi[(k+1)*dim+j]; 00441 00442 for (; dw > pts_per_bin; i++) { 00443 dw -= pts_per_bin; 00444 (xin[i]) = xnew - (xnew - xold) * dw / weight[k]; 00445 } 00446 } 00447 00448 for (k = 1 ; k < bins ; k++) { 00449 xi[k*dim+j] = xin[k]; 00450 } 00451 00452 xi[bins*dim+j]=1; 00453 } 00454 } 00455 return; 00456 } 00457 00458 /// Desc 00459 void print_lim(std::ostream *outs, 00460 const vec_t &xl, const vec_t &xu, unsigned long ldim) { 00461 unsigned long j; 00462 00463 (*outs) << "The limits of integration are:\n" << std::endl; 00464 for (j = 0; j < ldim; ++j) { 00465 (*outs) << "xl[" << j << "]=" << xl[j] << " xu[" << j 00466 << "]=" << xu[j] << std::endl; 00467 } 00468 (*outs) << std::endl; 00469 00470 return; 00471 } 00472 00473 /// Desc 00474 void print_head(std::ostream *outs, 00475 unsigned long num_dim, unsigned long calls, 00476 unsigned int lit_num, unsigned int lbins, 00477 unsigned int lboxes) { 00478 (*outs) << "num_dim=" << num_dim << ", calls=" << calls 00479 << ", it_num=" << lit_num << ", max_it_num=" 00480 << iterations << ", verb=" << this->verbose << ", alph=" << alpha 00481 << ",\n mode=" << mode << ", boxes=" << lboxes 00482 << "\n\n single.......iteration " 00483 << "accumulated......results\n" 00484 << "iter integral sigma integral " 00485 << " sigma chi-sq/it" << std::endl; 00486 00487 return; 00488 } 00489 00490 /// Desc 00491 void print_res(std::ostream *outs, unsigned int itr, double res, 00492 double err, double cum_res, double cum_err, 00493 double chi_sq) { 00494 (*outs).precision(5); 00495 (*outs) << itr << " "; 00496 outs->setf(std::ios::showpos); 00497 (*outs) << res << " "; 00498 outs->unsetf(std::ios::showpos); 00499 (*outs) << err << " "; 00500 outs->setf(std::ios::showpos); 00501 (*outs) << cum_res << " "; 00502 outs->unsetf(std::ios::showpos); 00503 (*outs) << cum_err << " " << chi_sq << std::endl; 00504 (*outs).precision(6); 00505 return; 00506 } 00507 00508 /// Desc 00509 void print_dist(std::ostream *outs, unsigned long ldim) { 00510 unsigned long i, j; 00511 int p = this->verbose; 00512 if (p < 1) 00513 return; 00514 00515 for (j = 0; j < ldim; ++j) { 00516 (*outs) << "\n axis " << j << std::endl; 00517 (*outs) << " x g" << std::endl; 00518 for (i = 0; i < bins; i++) { 00519 (*outs) << "weight [" << (xi[(i)*dim + (j)]) << " , " 00520 << xi[(i+1)*dim+j] << "] = "; 00521 (*outs) << " " << (d[(i)*dim + (j)]) << std::endl; 00522 00523 } 00524 (*outs) << std::endl; 00525 } 00526 (*outs) << std::endl; 00527 return; 00528 } 00529 00530 /// Desc 00531 void print_grid(std::ostream *outs, unsigned long ldim) { 00532 unsigned long i, j; 00533 int p = this->verbose; 00534 if (p < 1) { 00535 return; 00536 } 00537 00538 for (j = 0; j < ldim; ++j) { 00539 (*outs) << "\n axis " << j << std::endl; 00540 (*outs) << " x " << std::endl; 00541 for (i = 0; i <= bins; i++) { 00542 (*outs) << (xi[(i)*dim + (j)]) << " "; 00543 if (i % 5 == 4) { 00544 (*outs) << std::endl; 00545 } 00546 } 00547 (*outs) << std::endl; 00548 } 00549 (*outs) << std::endl; 00550 return; 00551 } 00552 00553 rng_t rng; 00554 00555 alloc_t ao; 00556 alloc_vec_t x; 00557 00558 public: 00559 00560 int vegas_integrate(const vec_t &xl, const vec_t &xu, 00561 size_t calls, double *lresult, double *abserr, 00562 func_t &func, param_t &pa) { 00563 00564 double cum_int, cum_sig; 00565 size_t i, k, it; 00566 00567 for (i = 0; i < dim; i++) { 00568 if (xu[i] <= xl[i]) { 00569 GSL_ERROR ("xu must be greater than xl", GSL_EINVAL); 00570 } 00571 00572 if (xu[i] - xl[i] > GSL_DBL_MAX) { 00573 GSL_ERROR ("Range of integration is too large, please rescale", 00574 GSL_EINVAL); 00575 } 00576 } 00577 00578 if (stage == 0) { 00579 init_grid (xl, xu, dim); 00580 00581 if ( this->verbose >= 0) { 00582 print_lim (&std::cout, xl, xu, dim); 00583 } 00584 } 00585 00586 if ( stage <= 1) { 00587 wtd_int_sum = 0; 00588 sum_wgts = 0; 00589 chi_sum = 0; 00590 it_num = 1; 00591 samples = 0; 00592 } 00593 00594 if ( stage <= 2) { 00595 unsigned int lbins = bins_max; 00596 unsigned int lboxes = 1; 00597 00598 if ( mode != GSL_VEGAS2_MODE_IMPORTANCE_ONLY) { 00599 /* shooting for 2 calls/box */ 00600 00601 lboxes = ((unsigned int)(floor (pow (calls / 2.0, 1.0 / dim)))); 00602 mode = GSL_VEGAS2_MODE_IMPORTANCE; 00603 00604 if (2 * lboxes >= bins_max) { 00605 /* if bins/box < 2 */ 00606 int box_per_bin = GSL_MAX (lboxes / bins_max, 1); 00607 00608 lbins = GSL_MIN(lboxes / box_per_bin, bins_max); 00609 lboxes = box_per_bin * lbins; 00610 00611 mode = GSL_VEGAS2_MODE_STRATIFIED; 00612 } 00613 00614 } 00615 00616 std::cout << calls << " " << dim << std::endl; 00617 std::cout << lboxes << " " << mode << " " << lbins 00618 << std::endl; 00619 00620 { 00621 double tot_boxes = pow ((double) lboxes, (double) dim); 00622 calls_per_box = ((unsigned int)(GSL_MAX (calls / 00623 tot_boxes, 2))); 00624 calls = ((size_t)( calls_per_box * tot_boxes)); 00625 } 00626 00627 /* total volume of x-space/(avg num of calls/bin) */ 00628 jac = vol * pow ((double) lbins, (double) dim) / calls; 00629 00630 /// FIXME check this 00631 lboxes = boxes; 00632 00633 /* If the number of bins changes from the previous invocation, bins 00634 are expanded or contracted accordingly, while preserving bin 00635 density */ 00636 00637 /// FIXME check this 00638 if (lbins != bins) { 00639 resize_grid (bins); 00640 00641 if ( this->verbose > 1) { 00642 print_grid(&std::cout,dim); 00643 } 00644 } 00645 00646 if ( this->verbose >= 0) { 00647 print_head(&std::cout,dim,calls, it_num, 00648 lbins, lboxes); 00649 } 00650 } 00651 00652 it_start = it_num; 00653 00654 cum_int = 0.0; 00655 cum_sig = 0.0; 00656 00657 chisq = 0.0; 00658 00659 for (it = 0; it < iterations; it++) { 00660 double intgrl = 0.0, intgrl_sq = 0.0; 00661 double sig = 0.0; 00662 double wgt; 00663 size_t lcalls_per_box = calls_per_box; 00664 double jacbin = jac; 00665 00666 it_num = it_start + it; 00667 00668 reset_grid_values (); 00669 init_box_coord ( box); 00670 00671 do { 00672 double m = 0, q = 0; 00673 double f_sq_sum = 0.0; 00674 00675 for (k = 0; k < lcalls_per_box; k++) { 00676 double fval, bin_vol; 00677 00678 random_point(x,bin,&bin_vol,box,xl,xu,rng); 00679 00680 func(dim,x,fval,pa); 00681 fval*=jacbin*bin_vol; 00682 00683 if (!finite(fval)) { 00684 std::cout << "NOt finite: " << std::endl; 00685 std::cout << x << std::endl; 00686 exit(-1); 00687 } 00688 00689 00690 /* recurrence for mean and variance */ 00691 00692 { 00693 double dt = fval - m; 00694 m += dt / (k + 1.0); 00695 q += dt * dt * (k / (k + 1.0)); 00696 } 00697 00698 if ( mode != GSL_VEGAS2_MODE_STRATIFIED) { 00699 double f_sq = fval * fval; 00700 accumulate_distribution (bin, f_sq); 00701 } 00702 } 00703 00704 intgrl += m * lcalls_per_box; 00705 00706 f_sq_sum = q * lcalls_per_box ; 00707 00708 sig += f_sq_sum ; 00709 00710 if ( mode == GSL_VEGAS2_MODE_STRATIFIED) { 00711 accumulate_distribution (bin, f_sq_sum); 00712 } 00713 } 00714 while ( change_box_coord ( box)); 00715 00716 /* Compute final results for this iteration */ 00717 00718 sig = sig / (lcalls_per_box - 1.0) ; 00719 00720 if (sig > 0) 00721 { 00722 wgt = 1.0 / sig; 00723 } 00724 else if ( sum_wgts > 0) 00725 { 00726 wgt = sum_wgts / samples; 00727 } 00728 else 00729 { 00730 wgt = 0.0; 00731 } 00732 00733 intgrl_sq = intgrl * intgrl; 00734 00735 /// FIXEM check this 00736 *lresult = intgrl; 00737 sigma = sqrt(sig); 00738 00739 if (wgt > 0.0) { 00740 samples++ ; 00741 sum_wgts += wgt; 00742 wtd_int_sum += intgrl * wgt; 00743 chi_sum += intgrl_sq * wgt; 00744 00745 cum_int = wtd_int_sum / sum_wgts; 00746 cum_sig = sqrt (1 / sum_wgts); 00747 00748 if ( samples > 1) { 00749 chisq = ( chi_sum - wtd_int_sum * 00750 cum_int) / ( samples - 1.0); 00751 } 00752 } else { 00753 cum_int += (intgrl - cum_int) / (it + 1.0); 00754 cum_sig = 0.0; 00755 } 00756 00757 00758 if ( this->verbose >= 0) { 00759 print_res (&std::cout, it_num, intgrl, sqrt (sig), 00760 cum_int, cum_sig, chisq); 00761 if (it + 1 == iterations && this->verbose > 0) { 00762 print_grid (&std::cout,dim); 00763 } 00764 } 00765 00766 if ( this->verbose > 1) { 00767 print_dist (&std::cout,dim); 00768 } 00769 00770 refine_grid (); 00771 00772 if ( this->verbose > 1) { 00773 print_grid (&std::cout,dim); 00774 } 00775 00776 } 00777 00778 /* By setting stage to 1 further calls will generate independent 00779 estimates based on the same grid, although it may be rebinned. */ 00780 00781 stage = 1; 00782 00783 *lresult = cum_int; 00784 *abserr = cum_sig; 00785 00786 return GSL_SUCCESS; 00787 } 00788 00789 virtual ~gsl_vegas2() {} 00790 00791 /// Integrate function \c func from x=a to x=b. 00792 virtual int minteg_err(func_t &func, size_t ndim, 00793 const vec_t &a, const vec_t &b, param_t &pa, 00794 double &res, double &err) 00795 { 00796 00797 ao.allocate(x,ndim); 00798 allocate(ndim); 00799 00800 int ncalls=this->n_points; 00801 00802 vegas_integrate(a,b,ncalls,&res,&err,func,pa); 00803 00804 free(); 00805 ao.free(x); 00806 00807 return 0; 00808 00809 } 00810 00811 /// Return string denoting type ("gsl_vegas2") 00812 virtual const char *type() { return "gsl_vegas2"; } 00813 00814 }; 00815 00816 #ifndef DOXYGENP 00817 } 00818 #endif 00819 00820 #endif 00821
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