1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767
|
@cindex Monte Carlo integration
@cindex stratified sampling in Monte Carlo integration
@cindex multidimensional integration
This chapter describes routines for multidimensional Monte Carlo
integration. These include the traditional Monte Carlo method and
adaptive algorithms such as @sc{vegas} and @sc{miser} which use
importance sampling and stratified sampling techniques. Each algorithm
computes an estimate of a multidimensional definite integral of the
form,
@tex
\beforedisplay
$$
I = \int_{x_l}^{x_u} dx\,\int_{y_l}^{y_u}dy\,... f(x,y,...)
$$
\afterdisplay
@end tex
@ifinfo
@example
I = \int_xl^xu dx \int_yl^yu dy ... f(x, y, ...)
@end example
@end ifinfo
@noindent
over a hypercubic region @math{((x_l,x_u)}, @math{(y_l,y_u), ...)} using
a fixed number of function calls. The routines also provide a
statistical estimate of the error on the result. This error estimate
should be taken as a guide rather than as a strict error bound---random
sampling of the region may not uncover all the important features
of the function, resulting in an underestimate of the error.
The functions are defined in separate header files for each routine,
@file{gsl_monte_plain.h}, @file{gsl_monte_miser.h} and
@file{gsl_monte_vegas.h}.
@menu
* Monte Carlo Interface::
* PLAIN Monte Carlo::
* MISER::
* VEGAS::
* Monte Carlo Examples::
* Monte Carlo Integration References and Further Reading::
@end menu
@node Monte Carlo Interface
@section Interface
All of the Monte Carlo integration routines use the same general form of
interface. There is an allocator to allocate memory for control
variables and workspace, a routine to initialize those control
variables, the integrator itself, and a function to free the space when
done.
Each integration function requires a random number generator to be
supplied, and returns an estimate of the integral and its standard
deviation. The accuracy of the result is determined by the number of
function calls specified by the user. If a known level of accuracy is
required this can be achieved by calling the integrator several times
and averaging the individual results until the desired accuracy is
obtained.
Random sample points used within the Monte Carlo routines are always
chosen strictly within the integration region, so that endpoint
singularities are automatically avoided.
The function to be integrated has its own datatype, defined in the
header file @file{gsl_monte.h}.
@deftp {Data Type} gsl_monte_function
This data type defines a general function with parameters for Monte
Carlo integration.
@table @code
@item double (* f) (double * @var{x}, size_t @var{dim}, void * @var{params})
this function should return the value
@c{$f(x,\hbox{\it params})$}
@math{f(x,params)} for the argument @var{x} and parameters @var{params},
where @var{x} is an array of size @var{dim} giving the coordinates of
the point where the function is to be evaluated.
@item size_t dim
the number of dimensions for @var{x}.
@item void * params
a pointer to the parameters of the function.
@end table
@end deftp
@noindent
Here is an example for a quadratic function in two dimensions,
@tex
\beforedisplay
$$
f(x,y) = a x^2 + b x y + c y^2
$$
\afterdisplay
@end tex
@ifinfo
@example
f(x,y) = a x^2 + b x y + c y^2
@end example
@end ifinfo
@noindent
with @math{a = 3}, @math{b = 2}, @math{c = 1}. The following code
defines a @code{gsl_monte_function} @code{F} which you could pass to an
integrator:
@example
struct my_f_params @{ double a; double b; double c; @};
double
my_f (double x[], size_t dim, void * p) @{
struct my_f_params * fp = (struct my_f_params *)p;
if (dim != 2)
@{
fprintf (stderr, "error: dim != 2");
abort ();
@}
return fp->a * x[0] * x[0]
+ fp->b * x[0] * x[1]
+ fp->c * x[1] * x[1];
@}
gsl_monte_function F;
struct my_f_params params = @{ 3.0, 2.0, 1.0 @};
F.f = &my_f;
F.dim = 2;
F.params = ¶ms;
@end example
@noindent
The function @math{f(x)} can be evaluated using the following macro,
@example
#define GSL_MONTE_FN_EVAL(F,x)
(*((F)->f))(x,(F)->dim,(F)->params)
@end example
@node PLAIN Monte Carlo
@section PLAIN Monte Carlo
@cindex plain Monte Carlo
The plain Monte Carlo algorithm samples points randomly from the
integration region to estimate the integral and its error. Using this
algorithm the estimate of the integral @math{E(f; N)} for @math{N}
randomly distributed points @math{x_i} is given by,
@tex
\beforedisplay
$$
E(f; N) = V \langle f \rangle = {V \over N} \sum_i^N f(x_i)
$$
\afterdisplay
@end tex
@ifinfo
@example
E(f; N) = = V <f> = (V / N) \sum_i^N f(x_i)
@end example
@end ifinfo
@noindent
where @math{V} is the volume of the integration region. The error on
this estimate @math{\sigma(E;N)} is calculated from the estimated
variance of the mean,
@tex
\beforedisplay
$$
\sigma^2 (E; N) = {V^2 \over N^2 } \sum_i^N (f(x_i) - \langle f \rangle)^2.
$$
\afterdisplay
@end tex
@ifinfo
@example
\sigma^2 (E; N) = (V^2 / N^2) \sum_i^N (f(x_i) - <f>)^2.
@end example
@end ifinfo
@noindent
For large @math{N} this variance decreases asymptotically as
@math{\Var(f)/N}, where @math{\Var(f)} is the true variance of the
function over the integration region. The error estimate itself should
decrease as @c{$\sigma(f)/\sqrt{N}$}
@math{\sigma(f)/\sqrt@{N@}}. The familiar law of errors
decreasing as @c{$1/\sqrt{N}$}
@math{1/\sqrt@{N@}} applies---to reduce the error by a
factor of 10 requires a 100-fold increase in the number of sample
points.
The functions described in this section are declared in the header file
@file{gsl_monte_plain.h}.
@deftypefun {gsl_monte_plain_state *} gsl_monte_plain_alloc (size_t @var{dim})
@tindex gsl_monte_plain_state
This function allocates and initializes a workspace for Monte Carlo
integration in @var{dim} dimensions.
@end deftypefun
@deftypefun int gsl_monte_plain_init (gsl_monte_plain_state* @var{s})
This function initializes a previously allocated integration state.
This allows an existing workspace to be reused for different
integrations.
@end deftypefun
@deftypefun int gsl_monte_plain_integrate (gsl_monte_function * @var{f}, const double @var{xl}[], const double @var{xu}[], size_t @var{dim}, size_t @var{calls}, gsl_rng * @var{r}, gsl_monte_plain_state * @var{s}, double * @var{result}, double * @var{abserr})
This routines uses the plain Monte Carlo algorithm to integrate the
function @var{f} over the @var{dim}-dimensional hypercubic region
defined by the lower and upper limits in the arrays @var{xl} and
@var{xu}, each of size @var{dim}. The integration uses a fixed number
of function calls @var{calls}, and obtains random sampling points using
the random number generator @var{r}. A previously allocated workspace
@var{s} must be supplied. The result of the integration is returned in
@var{result}, with an estimated absolute error @var{abserr}.
@end deftypefun
@deftypefun void gsl_monte_plain_free (gsl_monte_plain_state * @var{s})
This function frees the memory associated with the integrator state
@var{s}.
@end deftypefun
@node MISER
@section MISER
@cindex MISER monte carlo integration
@cindex recursive stratified sampling, MISER
The @sc{miser} algorithm of Press and Farrar is based on recursive
stratified sampling. This technique aims to reduce the overall
integration error by concentrating integration points in the regions of
highest variance.
The idea of stratified sampling begins with the observation that for two
disjoint regions @math{a} and @math{b} with Monte Carlo estimates of the
integral @math{E_a(f)} and @math{E_b(f)} and variances
@math{\sigma_a^2(f)} and @math{\sigma_b^2(f)}, the variance
@math{\Var(f)} of the combined estimate
@c{$E(f) = {1\over 2} (E_a(f) + E_b(f))$}
@math{E(f) = (1/2) (E_a(f) + E_b(f))}
is given by,
@tex
\beforedisplay
$$
\Var(f) = {\sigma_a^2(f) \over 4 N_a} + {\sigma_b^2(f) \over 4 N_b}.
$$
\afterdisplay
@end tex
@ifinfo
@example
\Var(f) = (\sigma_a^2(f) / 4 N_a) + (\sigma_b^2(f) / 4 N_b).
@end example
@end ifinfo
@noindent
It can be shown that this variance is minimized by distributing the
points such that,
@tex
\beforedisplay
$$
{N_a \over N_a+N_b} = {\sigma_a \over \sigma_a + \sigma_b}.
$$
\afterdisplay
@end tex
@ifinfo
@example
N_a / (N_a + N_b) = \sigma_a / (\sigma_a + \sigma_b).
@end example
@end ifinfo
@noindent
Hence the smallest error estimate is obtained by allocating sample
points in proportion to the standard deviation of the function in each
sub-region.
The @sc{miser} algorithm proceeds by bisecting the integration region
along one coordinate axis to give two sub-regions at each step. The
direction is chosen by examining all @math{d} possible bisections and
selecting the one which will minimize the combined variance of the two
sub-regions. The variance in the sub-regions is estimated by sampling
with a fraction of the total number of points available to the current
step. The same procedure is then repeated recursively for each of the
two half-spaces from the best bisection. The remaining sample points are
allocated to the sub-regions using the formula for @math{N_a} and
@math{N_b}. This recursive allocation of integration points continues
down to a user-specified depth where each sub-region is integrated using
a plain Monte Carlo estimate. These individual values and their error
estimates are then combined upwards to give an overall result and an
estimate of its error.
The functions described in this section are declared in the header file
@file{gsl_monte_miser.h}.
@deftypefun {gsl_monte_miser_state *} gsl_monte_miser_alloc (size_t @var{dim})
@tindex gsl_monte_miser_state
This function allocates and initializes a workspace for Monte Carlo
integration in @var{dim} dimensions. The workspace is used to maintain
the state of the integration.
@end deftypefun
@deftypefun int gsl_monte_miser_init (gsl_monte_miser_state* @var{s})
This function initializes a previously allocated integration state.
This allows an existing workspace to be reused for different
integrations.
@end deftypefun
@deftypefun int gsl_monte_miser_integrate (gsl_monte_function * @var{f}, const double @var{xl}[], const double @var{xu}[], size_t @var{dim}, size_t @var{calls}, gsl_rng * @var{r}, gsl_monte_miser_state * @var{s}, double * @var{result}, double * @var{abserr})
This routines uses the @sc{miser} Monte Carlo algorithm to integrate the
function @var{f} over the @var{dim}-dimensional hypercubic region
defined by the lower and upper limits in the arrays @var{xl} and
@var{xu}, each of size @var{dim}. The integration uses a fixed number
of function calls @var{calls}, and obtains random sampling points using
the random number generator @var{r}. A previously allocated workspace
@var{s} must be supplied. The result of the integration is returned in
@var{result}, with an estimated absolute error @var{abserr}.
@end deftypefun
@deftypefun void gsl_monte_miser_free (gsl_monte_miser_state * @var{s})
This function frees the memory associated with the integrator state
@var{s}.
@end deftypefun
The @sc{miser} algorithm has several configurable parameters which can
be changed using the following two functions.@footnote{The previous
method of accessing these fields directly through the
@code{gsl_monte_miser_state} struct is now deprecated.}
@deftypefun void gsl_monte_miser_params_get (const gsl_monte_miser_state * @var{s}, gsl_monte_miser_params * @var{params})
This function copies the parameters of the integrator state into the
user-supplied @var{params} structure.
@end deftypefun
@deftypefun void gsl_monte_miser_params_set (gsl_monte_miser_state * @var{s}, const gsl_monte_miser_params * @var{params})
This function sets the integrator parameters based on values provided
in the @var{params} structure.
@end deftypefun
Typically the values of the parameters are first read using
@code{gsl_monte_miser_params_get}, the necessary changes are made to
the fields of the @var{params} structure, and the values are copied
back into the integrator state using
@code{gsl_monte_miser_params_set}. The functions use the
@code{gsl_monte_miser_params} structure which contains the following
fields:
@deftypevar double estimate_frac
This parameter specifies the fraction of the currently available number of
function calls which are allocated to estimating the variance at each
recursive step. The default value is 0.1.
@end deftypevar
@deftypevar size_t min_calls
This parameter specifies the minimum number of function calls required
for each estimate of the variance. If the number of function calls
allocated to the estimate using @var{estimate_frac} falls below
@var{min_calls} then @var{min_calls} are used instead. This ensures
that each estimate maintains a reasonable level of accuracy. The
default value of @var{min_calls} is @code{16 * dim}.
@end deftypevar
@deftypevar size_t min_calls_per_bisection
This parameter specifies the minimum number of function calls required
to proceed with a bisection step. When a recursive step has fewer calls
available than @var{min_calls_per_bisection} it performs a plain Monte
Carlo estimate of the current sub-region and terminates its branch of
the recursion. The default value of this parameter is @code{32 *
min_calls}.
@end deftypevar
@deftypevar double alpha
This parameter controls how the estimated variances for the two
sub-regions of a bisection are combined when allocating points. With
recursive sampling the overall variance should scale better than
@math{1/N}, since the values from the sub-regions will be obtained using
a procedure which explicitly minimizes their variance. To accommodate
this behavior the @sc{miser} algorithm allows the total variance to
depend on a scaling parameter @math{\alpha},
@tex
\beforedisplay
$$
\Var(f) = {\sigma_a \over N_a^\alpha} + {\sigma_b \over N_b^\alpha}.
$$
\afterdisplay
@end tex
@ifinfo
@example
\Var(f) = @{\sigma_a \over N_a^\alpha@} + @{\sigma_b \over N_b^\alpha@}.
@end example
@end ifinfo
@noindent
The authors of the original paper describing @sc{miser} recommend the
value @math{\alpha = 2} as a good choice, obtained from numerical
experiments, and this is used as the default value in this
implementation.
@end deftypevar
@deftypevar double dither
This parameter introduces a random fractional variation of size
@var{dither} into each bisection, which can be used to break the
symmetry of integrands which are concentrated near the exact center of
the hypercubic integration region. The default value of dither is zero,
so no variation is introduced. If needed, a typical value of
@var{dither} is 0.1.
@end deftypevar
@node VEGAS
@section VEGAS
@cindex VEGAS Monte Carlo integration
@cindex importance sampling, VEGAS
The @sc{vegas} algorithm of Lepage is based on importance sampling. It
samples points from the probability distribution described by the
function @math{|f|}, so that the points are concentrated in the regions
that make the largest contribution to the integral.
In general, if the Monte Carlo integral of @math{f} is sampled with
points distributed according to a probability distribution described by
the function @math{g}, we obtain an estimate @math{E_g(f; N)},
@tex
\beforedisplay
$$
E_g(f; N) = E(f/g; N)
$$
\afterdisplay
@end tex
@ifinfo
@example
E_g(f; N) = E(f/g; N)
@end example
@end ifinfo
@noindent
with a corresponding variance,
@tex
\beforedisplay
$$
\Var_g(f; N) = \Var(f/g; N).
$$
\afterdisplay
@end tex
@ifinfo
@example
\Var_g(f; N) = \Var(f/g; N).
@end example
@end ifinfo
@noindent
If the probability distribution is chosen as @math{g = |f|/I(|f|)} then
it can be shown that the variance @math{V_g(f; N)} vanishes, and the
error in the estimate will be zero. In practice it is not possible to
sample from the exact distribution @math{g} for an arbitrary function, so
importance sampling algorithms aim to produce efficient approximations
to the desired distribution.
The @sc{vegas} algorithm approximates the exact distribution by making a
number of passes over the integration region while histogramming the
function @math{f}. Each histogram is used to define a sampling
distribution for the next pass. Asymptotically this procedure converges
to the desired distribution. In order
to avoid the number of histogram bins growing like @math{K^d} the
probability distribution is approximated by a separable function:
@c{$g(x_1, x_2,\ldots) = g_1(x_1) g_2(x_2)\ldots$}
@math{g(x_1, x_2, ...) = g_1(x_1) g_2(x_2) ...}
so that the number of bins required is only @math{Kd}.
This is equivalent to locating the peaks of the function from the
projections of the integrand onto the coordinate axes. The efficiency
of @sc{vegas} depends on the validity of this assumption. It is most
efficient when the peaks of the integrand are well-localized. If an
integrand can be rewritten in a form which is approximately separable
this will increase the efficiency of integration with @sc{vegas}.
@sc{vegas} incorporates a number of additional features, and combines both
stratified sampling and importance sampling. The integration region is
divided into a number of ``boxes'', with each box getting a fixed
number of points (the goal is 2). Each box can then have a fractional
number of bins, but if the ratio of bins-per-box is less than two, Vegas switches to a
kind variance reduction (rather than importance sampling).
@deftypefun {gsl_monte_vegas_state *} gsl_monte_vegas_alloc (size_t @var{dim})
@tindex gsl_monte_vegas_state
This function allocates and initializes a workspace for Monte Carlo
integration in @var{dim} dimensions. The workspace is used to maintain
the state of the integration.
@end deftypefun
@deftypefun int gsl_monte_vegas_init (gsl_monte_vegas_state* @var{s})
This function initializes a previously allocated integration state.
This allows an existing workspace to be reused for different
integrations.
@end deftypefun
@deftypefun int gsl_monte_vegas_integrate (gsl_monte_function * @var{f}, double @var{xl}[], double @var{xu}[], size_t @var{dim}, size_t @var{calls}, gsl_rng * @var{r}, gsl_monte_vegas_state * @var{s}, double * @var{result}, double * @var{abserr})
This routines uses the @sc{vegas} Monte Carlo algorithm to integrate the
function @var{f} over the @var{dim}-dimensional hypercubic region
defined by the lower and upper limits in the arrays @var{xl} and
@var{xu}, each of size @var{dim}. The integration uses a fixed number
of function calls @var{calls}, and obtains random sampling points using
the random number generator @var{r}. A previously allocated workspace
@var{s} must be supplied. The result of the integration is returned in
@var{result}, with an estimated absolute error @var{abserr}. The result
and its error estimate are based on a weighted average of independent
samples. The chi-squared per degree of freedom for the weighted average
is returned via the state struct component, @var{s->chisq}, and must be
consistent with 1 for the weighted average to be reliable.
@end deftypefun
@deftypefun void gsl_monte_vegas_free (gsl_monte_vegas_state * @var{s})
This function frees the memory associated with the integrator state
@var{s}.
@end deftypefun
The @sc{vegas} algorithm computes a number of independent estimates of the
integral internally, according to the @code{iterations} parameter
described below, and returns their weighted average. Random sampling of
the integrand can occasionally produce an estimate where the error is
zero, particularly if the function is constant in some regions. An
estimate with zero error causes the weighted average to break down and
must be handled separately. In the original Fortran implementations of
@sc{vegas} the error estimate is made non-zero by substituting a small
value (typically @code{1e-30}). The implementation in GSL differs from
this and avoids the use of an arbitrary constant---it either assigns
the value a weight which is the average weight of the preceding
estimates or discards it according to the following procedure,
@table @asis
@item current estimate has zero error, weighted average has finite error
The current estimate is assigned a weight which is the average weight of
the preceding estimates.
@item current estimate has finite error, previous estimates had zero error
The previous estimates are discarded and the weighted averaging
procedure begins with the current estimate.
@item current estimate has zero error, previous estimates had zero error
The estimates are averaged using the arithmetic mean, but no error is computed.
@end table
The convergence of the algorithm can be tested using the overall
chi-squared value of the results, which is available from the
following function:
@deftypefun double gsl_monte_vegas_chisq (const gsl_monte_vegas_state * @var{s})
This function returns the chi-squared per degree of freedom for the
weighted estimate of the integral. The returned value should be close
to 1. A value which differs significantly from 1 indicates that the
values from different iterations are inconsistent. In this case the
weighted error will be under-estimated, and further iterations of the
algorithm are needed to obtain reliable results.
@end deftypefun
@deftypefun void gsl_monte_vegas_runval (const gsl_monte_vegas_state * @var{s}, double * @var{result}, double * @var{sigma})
This function returns the raw (unaveraged) values of the integral
@var{result} and its error @var{sigma} from the most recent iteration
of the algorithm.
@end deftypefun
The @sc{vegas} algorithm is highly configurable. Several parameters
can be changed using the following two functions.
@deftypefun void gsl_monte_vegas_params_get (const gsl_monte_vegas_state * @var{s}, gsl_monte_vegas_params * @var{params})
This function copies the parameters of the integrator state into the
user-supplied @var{params} structure.
@end deftypefun
@deftypefun void gsl_monte_vegas_params_set (gsl_monte_vegas_state * @var{s}, const gsl_monte_vegas_params * @var{params})
This function sets the integrator parameters based on values provided
in the @var{params} structure.
@end deftypefun
Typically the values of the parameters are first read using
@code{gsl_monte_vegas_params_get}, the necessary changes are made to
the fields of the @var{params} structure, and the values are copied
back into the integrator state using
@code{gsl_monte_vegas_params_set}. The functions use the
@code{gsl_monte_vegas_params} structure which contains the following
fields:
@deftypevar double alpha
The parameter @code{alpha} controls the stiffness of the rebinning
algorithm. It is typically set between one and two. A value of zero
prevents rebinning of the grid. The default value is 1.5.
@end deftypevar
@deftypevar size_t iterations
The number of iterations to perform for each call to the routine. The
default value is 5 iterations.
@end deftypevar
@deftypevar int stage
Setting this determines the @dfn{stage} of the calculation. Normally,
@code{stage = 0} which begins with a new uniform grid and empty weighted
average. Calling @sc{vegas} with @code{stage = 1} retains the grid from the
previous run but discards the weighted average, so that one can ``tune''
the grid using a relatively small number of points and then do a large
run with @code{stage = 1} on the optimized grid. Setting @code{stage =
2} keeps the grid and the weighted average from the previous run, but
may increase (or decrease) the number of histogram bins in the grid
depending on the number of calls available. Choosing @code{stage = 3}
enters at the main loop, so that nothing is changed, and is equivalent
to performing additional iterations in a previous call.
@end deftypevar
@deftypevar int mode
The possible choices are @code{GSL_VEGAS_MODE_IMPORTANCE},
@code{GSL_VEGAS_MODE_STRATIFIED}, @code{GSL_VEGAS_MODE_IMPORTANCE_ONLY}.
This determines whether @sc{vegas} will use importance sampling or
stratified sampling, or whether it can pick on its own. In low
dimensions @sc{vegas} uses strict stratified sampling (more precisely,
stratified sampling is chosen if there are fewer than 2 bins per box).
@end deftypevar
@deftypevar int verbose
@deftypevarx {FILE *} ostream
These parameters set the level of information printed by @sc{vegas}. All
information is written to the stream @var{ostream}. The default setting
of @var{verbose} is @code{-1}, which turns off all output. A
@var{verbose} value of @code{0} prints summary information about the
weighted average and final result, while a value of @code{1} also
displays the grid coordinates. A value of @code{2} prints information
from the rebinning procedure for each iteration.
@end deftypevar
The above fields and the @var{chisq} value can also be accessed
directly in the @code{gsl_monte_vegas_state} but such use is
deprecated.
@node Monte Carlo Examples
@section Examples
The example program below uses the Monte Carlo routines to estimate the
value of the following 3-dimensional integral from the theory of random
walks,
@tex
\beforedisplay
$$
I = \int_{-\pi}^{+\pi} {dk_x \over 2\pi}
\int_{-\pi}^{+\pi} {dk_y \over 2\pi}
\int_{-\pi}^{+\pi} {dk_z \over 2\pi}
{ 1 \over (1 - \cos(k_x)\cos(k_y)\cos(k_z))}.
$$
\afterdisplay
@end tex
@ifinfo
@example
I = \int_@{-pi@}^@{+pi@} @{dk_x/(2 pi)@}
\int_@{-pi@}^@{+pi@} @{dk_y/(2 pi)@}
\int_@{-pi@}^@{+pi@} @{dk_z/(2 pi)@}
1 / (1 - cos(k_x)cos(k_y)cos(k_z)).
@end example
@end ifinfo
@noindent
The analytic value of this integral can be shown to be @math{I =
\Gamma(1/4)^4/(4 \pi^3) = 1.393203929685676859...}. The integral gives
the mean time spent at the origin by a random walk on a body-centered
cubic lattice in three dimensions.
For simplicity we will compute the integral over the region
@math{(0,0,0)} to @math{(\pi,\pi,\pi)} and multiply by 8 to obtain the
full result. The integral is slowly varying in the middle of the region
but has integrable singularities at the corners @math{(0,0,0)},
@math{(0,\pi,\pi)}, @math{(\pi,0,\pi)} and @math{(\pi,\pi,0)}. The
Monte Carlo routines only select points which are strictly within the
integration region and so no special measures are needed to avoid these
singularities.
@smallexample
@verbatiminclude examples/monte.c
@end smallexample
@noindent
With 500,000 function calls the plain Monte Carlo algorithm achieves a
fractional error of 1%. The estimated error @code{sigma} is roughly
consistent with the actual error--the computed result differs from
the true result by about 1.4 standard deviations,
@example
plain ==================
result = 1.412209
sigma = 0.013436
exact = 1.393204
error = 0.019005 = 1.4 sigma
@end example
@noindent
The @sc{miser} algorithm reduces the error by a factor of four, and also
correctly estimates the error,
@example
miser ==================
result = 1.391322
sigma = 0.003461
exact = 1.393204
error = -0.001882 = 0.54 sigma
@end example
@noindent
In the case of the @sc{vegas} algorithm the program uses an initial
warm-up run of 10,000 function calls to prepare, or ``warm up'', the grid.
This is followed by a main run with five iterations of 100,000 function
calls. The chi-squared per degree of freedom for the five iterations are
checked for consistency with 1, and the run is repeated if the results
have not converged. In this case the estimates are consistent on the
first pass.
@example
vegas warm-up ==================
result = 1.392673
sigma = 0.003410
exact = 1.393204
error = -0.000531 = 0.16 sigma
converging...
result = 1.393281 sigma = 0.000362 chisq/dof = 1.5
vegas final ==================
result = 1.393281
sigma = 0.000362
exact = 1.393204
error = 0.000077 = 0.21 sigma
@end example
@noindent
If the value of @code{chisq} had differed significantly from 1 it would
indicate inconsistent results, with a correspondingly underestimated
error. The final estimate from @sc{vegas} (using a similar number of
function calls) is significantly more accurate than the other two
algorithms.
@node Monte Carlo Integration References and Further Reading
@section References and Further Reading
The @sc{miser} algorithm is described in the following article by Press
and Farrar,
@itemize @w{}
@item
W.H. Press, G.R. Farrar, @cite{Recursive Stratified Sampling for
Multidimensional Monte Carlo Integration},
Computers in Physics, v4 (1990), pp190--195.
@end itemize
@noindent
The @sc{vegas} algorithm is described in the following papers,
@itemize @w{}
@item
G.P. Lepage,
@cite{A New Algorithm for Adaptive Multidimensional Integration},
Journal of Computational Physics 27, 192--203, (1978)
@item
G.P. Lepage,
@cite{VEGAS: An Adaptive Multi-dimensional Integration Program},
Cornell preprint CLNS 80-447, March 1980
@end itemize
|