vegas Module

Introduction

The key Python objects supported by the vegas module are:

  • vegas.Integrator — an object describing a multidimensional integration operator. Such objects contain information about the integration volume, and also about optimal remappings of the integration variables based upon the last integral evaluated using the object.

  • vegas.AdaptiveMap — an object describing the remappings used by vegas.

  • vegas.RAvg — an object describing the result of a vegas integration. vegas returns the weighted average of the integral estimates from each vegas iteration as an object of class vegas.RAvg. These are Gaussian random variables — that is, they have a mean and a standard deviation — but also contain information about the iterations vegas used in generating the result.

  • vegas.RAvgArray — an array version of vegas.RAvg used when the integrand is array-valued.

  • vegas.RAvgDict — a dictionary version of vegas.RAvg used when the integrand is dictionary-valued.

  • vegas.PDFIntegrator — a specialized integrator for evaluating Gaussian expectation values.

These are described in detail below.

Integrator Objects

The central component of the vegas package is the integrator class:

class vegas.Integrator

Adaptive multidimensional Monte Carlo integration.

vegas.Integrator objects make Monte Carlo estimates of multidimensional functions f(x) where x[d] is a point in the integration volume:

integ = vegas.Integrator(integration_region)

result = integ(f, nitn=10, neval=10000)

The integator makes nitn estimates of the integral, each using at most neval samples of the integrand, as it adapts to the specific features of the integrand. Successive estimates (iterations) typically improve in accuracy until the integrator has fully adapted. The integrator returns the weighted average of all nitn estimates, together with an estimate of the statistical (Monte Carlo) uncertainty in that estimate of the integral. The result is an object of type RAvg (which is derived from gvar.GVar).

Integrands f(x) return numbers, arrays of numbers (any shape), or dictionaries whose values are numbers or arrays (any shape). Each number returned by an integrand corresponds to a different integrand. When arrays are returned, vegas adapts to the first number in the flattened array. When dictionaries are returned, vegas adapts to the first number in the value corresponding to the first key.

vegas can generate integration points in batches for integrands built from classes derived from vegas.BatchIntegrand, or integrand functions decorated by vegas.batchintegrand(). Batch integrands are typically much faster, especially if they are coded in Cython.

vegas.Integrators have a large number of parameters but the only ones that most people will care about are: the number nitn of iterations of the vegas algorithm; the maximum number neval of integrand evaluations per iteration; and the damping parameter alpha, which is used to slow down the adaptive algorithms when they would otherwise be unstable (e.g., with very peaky integrands). Setting parameter analyzer=vegas.reporter() is sometimes useful, as well, since it causes vegas to print (on sys.stdout) intermediate results from each iteration, as they are produced. This helps when each iteration takes a long time to complete (e.g., longer than an hour) because it allows you to monitor progress as it is being made (or not).

Parameters
  • map (array, vegas.AdaptiveMap or vegas.Integrator) –

    The integration region as specified by an array map[d, i] where d is the direction and i=0,1 specify the lower and upper limits of integration in direction d.

    map could also be the integration map from another vegas.Integrator, or that vegas.Integrator itself. In this case the grid is copied from the existing integrator.

  • nitn (positive int) – The maximum number of iterations used to adapt to the integrand and estimate its value. The default value is 10; typical values range from 10 to 20.

  • neval (positive int) – Approximate number of integrand evaluations in each iteration of the vegas algorithm. Increasing neval increases the precision: statistical errors should fall at least as fast as sqrt(1./neval) and often fall much faster. When beta>0, the total number of evaluations can be substantially larger than neval for challenging integrands; parameter max_neval_hcube can be used to limit such growth. The default value is 1000; real problems often require 10–1000 times more evaluations than this.

  • alpha (float) – Damping parameter controlling the remapping of the integration variables as vegas adapts to the integrand. Smaller values slow adaptation, which may be desirable for difficult integrands. Small or zero alphas are also sometimes useful after the grid has adapted, to minimize fluctuations away from the optimal grid. The default value is 0.5.

  • beta (float) – Damping parameter controlling the redistribution of integrand evaluations across hypercubes in the stratified sampling of the integral (over transformed variables). Smaller values limit the amount of redistribution. The theoretically optimal value is 1; setting beta=0 prevents any redistribution of evaluations. The default value is 0.75.

  • adapt (bool) – Setting adapt=False prevents further adaptation by vegas. Typically this would be done after training the vegas.Integrator on an integrand, in order to stabilize further estimates of the integral. vegas uses unweighted averages to combine results from different iterations when adapt=False. The default setting is adapt=True.

  • nhcube_batch (positive int) – The number of hypercubes (in y space) whose integration points are combined into a single batch to be passed to the integrand, together, when using vegas in batch mode. The default value is 1000. Larger values may be lead to faster evaluations, but at the cost of more memory for internal work arrays.

  • minimize_mem (bool) – When True, vegas minimizes internal workspace at the cost of extra evaluations of the integrand. This can increase execution time by 50–100% but might be desirable when the number of evaluations is very large (e.g., neval=1e9). Normally vegas uses internal work space that grows in proportion to neval. If that work space exceeds the size of the RAM available to the processor, vegas runs much more slowly. Setting minimize_mem=True greatly reduces the internal storage used by vegas; in particular memory becomes independent of neval. The default setting (minimize_mem=False), however, is much superior unless memory becomes a problem. (The large memory is needed for adaptive stratified sampling, so memory is not an issue if beta=0.)

  • adapt_to_errors (bool) –

    adapt_to_errors=False causes vegas to remap the integration variables to emphasize regions where |f(x)| is largest. This is the default mode.

    adapt_to_errors=True causes vegas to remap variables to emphasize regions where the Monte Carlo error is largest. This might be superior when the number of the number of stratifications (self.nstrat) in the y grid is large (> 50?). It is typically useful only in one or two dimensions.

  • maxinc_axis (positive int) – The maximum number of increments per axis allowed for the x-space grid. The default value is 1000; there is probably little need to use other values.

  • max_nhcube (positive int) – Maximum number of y-space hypercubes used for stratified sampling. Setting max_nhcube=1 turns stratified sampling off, which is probably never a good idea. The default setting (1e9) was chosen to correspond to the point where internal work arrays become comparable in size to the typical amount of RAM available to a processor (in a laptop in 2014). Internal memory usage is large only when beta>0 and minimize_mem=False; therefore max_nhcube is ignored if beta=0 or minimize_mem=True.

  • max_neval_hcube (positive int) – Maximum number of integrand evaluations per hypercube in the stratification. The default value is 1e7. Larger values might allow for more adaptation (when beta>0), but also can result in very large internal work arrays.

  • rtol (float) – Relative error in the integral estimate at which point the integrator can stop. The default value is 0.0 which turns off this stopping condition. This stopping condition can be quite unreliable in early iterations, before vegas has converged. Use with caution, if at all.

  • atol (float) – Absolute error in the integral estimate at which point the integrator can stop. The default value is 0.0 which turns off this stopping condition. This stopping condition can be quite unreliable in early iterations, before vegas has converged. Use with caution, if at all.

  • analyzer

    An object with methods

    analyzer.begin(itn, integrator)

    analyzer.end(itn_result, result)

    where: begin(itn, integrator) is called at the start of each vegas iteration with itn equal to the iteration number and integrator equal to the integrator itself; and end(itn_result, result) is called at the end of each iteration with itn_result equal to the result for that iteration and result equal to the cummulative result of all iterations so far. Setting analyzer=vegas.reporter(), for example, causes vegas to print out a running report of its results as they are produced. The default is analyzer=None.

  • ran_array_generator – Function that generates numpy arrays of random numbers distributed uniformly between 0 and 1. ran_array_generator(shape) should create an array whose dimensions are specified by the integer-valued tuple shape. The default generator is numpy.random.random.

  • sync_ran (bool) – If True, the default random number generator is synchronized across all processors when using MPI. If False, vegas does no synchronization (but the random numbers should synchronized some other way).

  • mpi – Setting mpi=False disables mpi support in vegas even if mpi is available; setting mpi=True (default) allows use of mpi provided module mpi4py is installed.

Methods:
__call__(fcn, **kargs)

Integrate integrand fcn.

A typical integrand has the form, for example:

def f(x):
    return x[0] ** 2 + x[1] ** 4

The argument x[d] is an integration point, where index d=0... represents direction within the integration volume.

Integrands can be array-valued, representing multiple integrands: e.g.,

def f(x):
    return [x[0] ** 2, x[0] / x[1]]

The return arrays can have any shape. Dictionary-valued integrands are also supported: e.g.,

def f(x):
    return dict(a=x[0] ** 2, b=[x[0] / x[1], x[1] / x[0]])

Integrand functions that return arrays or dictionaries are useful for multiple integrands that are closely related, and can lead to substantial reductions in the errors for ratios or differences of the results.

It is usually much faster to use vegas in batch mode, where integration points are presented to the integrand in batches. A simple batch integrand might be, for example:

@vegas.batchintegrand
def f(x):
    return x[:, 0] ** 2 + x[:, 1] ** 4

where decorator @vegas.batchintegrand tells vegas that the integrand processes integration points in batches. The array x[i, d] represents a collection of different integration points labeled by i=0.... (The number of points is controlled vegas.Integrator parameter nhcube_batch.) The batch index is always first.

Batch integrands can also be constructed from classes derived from vegas.BatchIntegrand.

Batch mode is particularly useful (and fast) when the class derived from vegas.BatchIntegrand is coded in Cython. Then loops over the integration points can be coded explicitly, avoiding the need to use numpy’s whole-array operators if they are not well suited to the integrand.

Any vegas parameter can also be reset: e.g., self(fcn, nitn=20, neval=1e6).

set(ka={}, **kargs)

Reset default parameters in integrator.

Usage is analogous to the constructor for vegas.Integrator: for example,

old_defaults = integ.set(neval=1e6, nitn=20)

resets the default values for neval and nitn in vegas.Integrator integ. A dictionary, here old_defaults, is returned. It can be used to restore the old defaults using, for example:

integ.set(old_defaults)
settings(ngrid=0)

Assemble summary of integrator settings into string.

Parameters

ngrid (int) – Number of grid nodes in each direction to include in summary. The default is 0.

Returns

String containing the settings.

random(yield_hcube=False, yield_y=False)

Iterator over integration points and weights.

This method creates an iterator that returns integration points from vegas, and their corresponding weights in an integral. Each point x[d] is accompanied by the weight assigned to that point by vegas when estimating an integral. Optionally it will also return the index of the hypercube containing the integration point and/or the y-space coordinates:

integ.random()  yields  x, wgt

integ.random(yield_hcube=True) yields x, wgt, hcube

integ.random(yield_y=True) yields x, y, wgt

integ.random(yield_hcube=True, yield_y=True) yields x, y, wgt, hcube

The number of integration points returned by the iterator corresponds to a single iteration.

random_batch(yield_hcube=False, yield_y=False)

Iterator over integration points and weights.

This method creates an iterator that returns integration points from vegas, and their corresponding weights in an integral. The points are provided in arrays x[i, d] where i=0... labels the integration points in a batch and d=0... labels direction. The corresponding weights assigned by vegas to each point are provided in an array wgt[i].

Optionally the integrator will also return the indices of the hypercubes containing the integration points and/or the y-space coordinates of those points:

integ.random_batch()  yields  x, wgt

integ.random_batch(yield_hcube=True) yields x, wgt, hcube

integ.random_batch(yield_y=True) yields x, y, wgt

integ.random_batch(yield_hcube=True, yield_y=True) yields x, y, wgt, hcube

The number of integration points returned by the iterator corresponds to a single iteration. The number in a batch is controlled by parameter nhcube_batch.

AdaptiveMap Objects

vegas’s remapping of the integration variables is handled by a vegas.AdaptiveMap object, which maps the original integration variables x into new variables y in a unit hypercube. Each direction has its own map specified by a grid in x space:

x_0 &= a \\
x_1 &= x_0 + \Delta x_0 \\
x_2 &= x_1 + \Delta x_1 \\
\cdots \\
x_N &= x_{N-1} + \Delta x_{N-1} = b

where a and b are the limits of integration. The grid specifies the transformation function at the points y=i/N for i=0,1\ldots N:

x(y\!=\!i/N) = x_i

Linear interpolation is used between those points. The Jacobian for this transformation is:

J(y) = J_i = N \Delta x_i

vegas adjusts the increments sizes to optimize its Monte Carlo estimates of the integral. This involves training the grid. To illustrate how this is done with vegas.AdaptiveMaps consider a simple two dimensional integral over a unit hypercube with integrand:

def f(x):
   return x[0] * x[1] ** 2

We want to create a grid that optimizes uniform Monte Carlo estimates of the integral in y space. We do this by sampling the integrand at a large number ny of random points y[j, d], where j=0...ny-1 and d=0,1, uniformly distributed throughout the integration volume in y space. These samples be used to train the grid using the following code:

import vegas
import numpy as np

def f(x):
   return x[0] * x[1] ** 2

m = vegas.AdaptiveMap([[0, 1], [0, 1]], ninc=5)

ny = 1000
y = np.random.uniform(0., 1., (ny, 2))  # 1000 random y's

x = np.empty(y.shape, float)            # work space
jac = np.empty(y.shape[0], float)
f2 = np.empty(y.shape[0], float)

print('intial grid:')
print(m.settings())

for itn in range(5):                    # 5 iterations to adapt
   m.map(y, x, jac)                     # compute x's and jac

   for j in range(ny):                  # compute training data
      f2[j] = (jac[j] * f(x[j])) ** 2

   m.add_training_data(y, f2)           # adapt
   m.adapt(alpha=1.5)

   print('iteration %d:' % itn)
   print(m.settings())

In each of the 5 iterations, the vegas.AdaptiveMap adjusts the map, making increments smaller where f2 is larger and larger where f2 is smaller. The map converges after only 2 or 3 iterations, as is clear from the output:

initial grid:
    grid[ 0] = [ 0.   0.2  0.4  0.6  0.8  1. ]
    grid[ 1] = [ 0.   0.2  0.4  0.6  0.8  1. ]

iteration 0:
    grid[ 0] = [ 0.     0.412  0.62   0.76   0.883  1.   ]
    grid[ 1] = [ 0.     0.506  0.691  0.821  0.91   1.   ]

iteration 1:
    grid[ 0] = [ 0.     0.428  0.63   0.772  0.893  1.   ]
    grid[ 1] = [ 0.     0.531  0.713  0.832  0.921  1.   ]

iteration 2:
    grid[ 0] = [ 0.     0.433  0.63   0.772  0.894  1.   ]
    grid[ 1] = [ 0.     0.533  0.714  0.831  0.922  1.   ]

iteration 3:
    grid[ 0] = [ 0.     0.435  0.631  0.772  0.894  1.   ]
    grid[ 1] = [ 0.     0.533  0.715  0.831  0.923  1.   ]

iteration 4:
    grid[ 0] = [ 0.     0.436  0.631  0.772  0.895  1.   ]
    grid[ 1] = [ 0.     0.533  0.715  0.831  0.924  1.   ]

The grid increments along direction 0 shrink at larger values x[0], varying as 1/x[0]. Along direction 1 the increments shrink more quickly varying like 1/x[1]**2.

vegas samples the integrand in order to estimate the integral. It uses those same samples to train its vegas.AdaptiveMap in this fashion, for use in subsequent iterations of the algorithm.

class vegas.AdaptiveMap

Adaptive map y->x(y) for multidimensional y and x.

An AdaptiveMap defines a multidimensional map y -> x(y) from the unit hypercube, with 0 <= y[d] <= 1, to an arbitrary hypercube in x space. Each direction is mapped independently with a Jacobian that is tunable (i.e., “adaptive”).

The map is specified by a grid in x-space that, by definition, maps into a uniformly spaced grid in y-space. The nodes of the grid are specified by grid[d, i] where d is the direction (d=0,1...dim-1) and i labels the grid point (i=0,1...N). The mapping for a specific point y into x space is:

y[d] -> x[d] = grid[d, i(y[d])] + inc[d, i(y[d])] * delta(y[d])

where i(y)=floor(y*N), delta(y)=y*N - i(y), and inc[d, i] = grid[d, i+1] - grid[d, i]. The Jacobian for this map,

dx[d]/dy[d] = inc[d, i(y[d])] * N,

is piece-wise constant and proportional to the x-space grid spacing. Each increment in the x-space grid maps into an increment of size 1/N in the corresponding y space. So regions in x space where inc[d, i] is small are stretched out in y space, while larger increments are compressed.

The x grid for an AdaptiveMap can be specified explicitly when the map is created: for example,

m = AdaptiveMap([[0, 0.1, 1], [-1, 0, 1]])

creates a two-dimensional map where the x[0] interval (0,0.1) and (0.1,1) map into the y[0] intervals (0,0.5) and (0.5,1) respectively, while x[1] intervals (-1,0) and (0,1) map into y[1] intervals (0,0.5) and (0.5,1).

More typically an initially uniform map is trained with data f[j] corresponding to ny points y[j, d], with j=0...ny-1, uniformly distributed in y space: for example,

m.add_training_data(y, f)
m.adapt(alpha=1.5)

m.adapt(alpha=1.5) shrinks grid increments where f[j] is large, and expands them where f[j] is small. Typically one has to iterate over several sets of ys and fs before the grid has fully adapted.

The speed with which the grid adapts is determined by parameter alpha. Large (positive) values imply rapid adaptation, while small values (much less than one) imply slow adaptation. As in any iterative process, it is usually a good idea to slow adaptation down in order to avoid instabilities.

Parameters
  • grid (array) – Initial x grid, where grid[d, i] is the i-th node in direction d.

  • ninc (int or None) – Number of increments along each axis of the x grid. A new grid is generated if ninc differs from grid.shape[1]. The new grid is designed to give the same Jacobian dx(y)/dy as the original grid. The default value, ninc=None, leaves the grid unchanged.

Attributes and methods:
dim

Number of dimensions.

ninc

Number of increments along each grid axis.

grid

The nodes of the grid defining the maps are self.grid[d, i] where d=0... specifies the direction and i=0...self.ninc the node.

inc

The increment widths of the grid:

self.inc[d, i] = self.grid[d, i + 1] - self.grid[d, i]
adapt(alpha=0.0, ninc=None)

Adapt grid to accumulated training data.

self.adapt(...) projects the training data onto each axis independently and maps it into x space. It shrinks x-grid increments in regions where the projected training data is large, and grows increments where the projected data is small. The grid along any direction is unchanged if the training data is constant along that direction.

The number of increments along a direction can be changed by setting parameter ninc.

The grid does not change if no training data has been accumulated, unless ninc is specified, in which case the number of increments is adjusted while preserving the relative density of increments at different values of x.

Parameters
  • alpha (double) – Determines the speed with which the grid adapts to training data. Large (postive) values imply rapid evolution; small values (much less than one) imply slow evolution. Typical values are of order one. Choosing alpha<0 causes adaptation to the unmodified training data (usually not a good idea).

  • ninc (int) – Number of increments along each direction in the new grid. The number is unchanged from the old grid if ninc is omitted (or equals None).

add_training_data(y, f, ny=-1)

Add training data f for y-space points y.

Accumulates training data for later use by self.adapt(). Grid increments will be made smaller in regions where f is larger than average, and larger where f is smaller than average. The grid is unchanged (converged?) when f is constant across the grid.

Parameters
  • y (array) – y values corresponding to the training data. y is a contiguous 2-d array, where y[j, d] is for points along direction d.

  • f (array) – Training function values. f[j] corresponds to point y[j, d] in y-space.

  • ny (int) – Number of y points: y[j, d] for d=0...dim-1 and j=0...ny-1. ny is set to y.shape[0] if it is omitted (or negative).

__call__(y)

Return x values corresponding to y.

y can be a single dim-dimensional point, or it can be an array y[i,j, ..., d] of such points (d=0..dim-1).

If y=None (default), y is set equal to a (uniform) random point in the volume.

jac(y)

Return the map’s Jacobian at y.

y can be a single dim-dimensional point, or it can be an array y[d,i,j,...] of such points (d=0..dim-1).

make_uniform(ninc=None)

Replace the grid with a uniform grid.

The new grid has ninc increments along each direction if ninc is specified. Otherwise it has the same number of increments as the old grid.

map(y, x, jac, ny=-1)

Map y to x, where jac is the Jacobian.

y[j, d] is an array of ny y-values for direction d. x[j, d] is filled with the corresponding x values, and jac[j] is filled with the corresponding Jacobian values. x and jac must be preallocated: for example,

x = numpy.empty(y.shape, numpy.float_)
jac = numpy.empty(y.shape[0], numpy.float_)
Parameters
  • y (array) – y values to be mapped. y is a contiguous 2-d array, where y[j, d] contains values for points along direction d.

  • x (array) – Container for x[j, d] values corresponding to y[j, d].

  • jac (array) – Container for Jacobian values jac[j] corresponding to ``y[j, d].

  • ny (int) – Number of y points: y[j, d] for d=0...dim-1 and j=0...ny-1. ny is set to y.shape[0] if it is omitted (or negative).

show_grid(ngrid=40, shrink=False)

Display plots showing the current grid.

Parameters
  • ngrid (int) – The number of grid nodes in each direction to include in the plot. The default is 40.

  • axes – List of pairs of directions to use in different views of the grid. Using None in place of a direction plots the grid for only one direction. Omitting axes causes a default set of pairings to be used.

  • shrink – Display entire range of each axis if False; otherwise shrink range to include just the nodes being displayed. The default is False.

settings(ngrid=5)

Create string with information about grid nodes.

Creates a string containing the locations of the nodes in the map grid for each direction. Parameter ngrid specifies the maximum number of nodes to print (spread evenly over the grid).

PDFIntegrator Objects

Expectation values using probability density functions defined by collections of Gaussian random variables (see gvar) can be evaluated using the following specialized integrator:

class vegas.PDFIntegrator(g, limit=1e15, scale=1., svdcut=1e-15)

vegas integrator for PDF expectation values.

PDFIntegrator(g) is a vegas integrator that evaluates expectation values for the multi-dimensional Gaussian distribution specified by with g, which is a gvar.GVar or an array of gvar.GVars or a dictionary whose values are gvar.GVars or arrays of gvar.GVars.

PDFIntegrator integrates over the entire parameter space of the distribution but reexpresses integrals in terms of variables that diagonalize g’s covariance matrix and are centered at its mean value. This greatly facilitates integration over these variables using the vegas module, making integrals over 10s or more of parameters feasible.

A simple illustration of PDFIntegrator is given by the following code:

import vegas
import gvar as gv

# multi-dimensional Gaussian distribution
g = gv.BufferDict()
g['a'] = gv.gvar([0., 1.], [[1., 0.9], [0.9, 1.]])
g['b'] = gv.gvar('1(1)')

# integrator for expectation values in distribution g
g_expval = vegas.PDFIntegrator(g)

# want expectation value of [fp, fp**2]
def f_f2(p):
    fp = p['a'][0] * p['a'][1] + p['b']
    return [fp, fp ** 2]

# adapt integrator to f_f2
warmup = g_expval(f_f2, neval=1000, nitn=5)

# <f_f2> in distribution g
results = g_expval(f_f2, neval=1000, nitn=5, adapt=False)
print(results.summary())
print('results =', results, '\n')

# mean and standard deviation of f(p)'s distribution
fmean = results[0]
fsdev = gv.sqrt(results[1] - results[0] ** 2)
print ('f.mean =', fmean, '   f.sdev =', fsdev)
print ("Gaussian approx'n for f(g) =", f_f2(g)[0])

where the warmup calls to the integrator are used to adapt it to the integrand, and the final results are in results. Here neval is the (approximate) number of function calls per iteration of the vegas algorithm and nitn is the number of iterations. We use the integrator to calculated the expectation value of fp and fp**2, so we can compute the standard deviation for the distribution of fps. The output from this code shows that the Gaussian approximation (1.0(1.4)) for the mean and standard deviation of the fp distribution is not particularly accurate here (correct value is 1.9(2.0)), because of the large uncertainties in g:

itn   integral        average         chi2/dof        Q
-------------------------------------------------------
  1   1.893(38)       1.893(38)           0.00     1.00
  2   1.905(35)       1.899(26)           0.25     0.78
  3   1.854(41)       1.884(22)           0.47     0.76
  4   1.921(36)       1.893(19)           0.44     0.85
  5   1.913(37)       1.897(17)           0.35     0.94

results = [1.897(17) 7.48(10)]

f.mean = 1.897(17)    f.sdev = 1.969(21)
Gaussian approx'n for f(g) = 1.0(1.4)

In general functions being integrated can return a number, or an array of numbers, or a dictionary whose values are numbers or arrays of numbers. This allows multiple expectation values to be evaluated simultaneously.

See the documentation with the vegas module for more details on its use, and on the attributes and methods associated with integrators. The example above sets adapt=False when computing final results. This gives more reliable error estimates when neval is small. Note that neval may need to be much larger (tens or hundreds of thousands) for more difficult high-dimension integrals.

Parameters
  • ggvar.GVar, array of gvar.GVars, or dictionary whose values are gvar.GVars or arrays of gvar.GVars that specifies the multi-dimensional Gaussian distribution used to construct the probability density function.

  • limit (positive float) – Limits the integrations to a finite region of size limit times the standard deviation on either side of the mean. This can be useful if the functions being integrated misbehave for large parameter values (e.g., numpy.exp overflows for a large range of arguments). Default is 1e15.

  • scale (positive float) – The integration variables are rescaled to emphasize parameter values of order scale times the standard deviation. The rescaling does not change the value of the integral but it can reduce uncertainties in the vegas estimate. Default is 1.0.

  • svdcut (non-negative float or None) – If not None, replace covariance matrix of g with a new matrix whose small eigenvalues are modified: eigenvalues smaller than svdcut times the maximum eigenvalue eig_max are replaced by svdcut*eig_max. This can ameliorate problems caused by roundoff errors when inverting the covariance matrix. It increases the uncertainty associated with the modified eigenvalues and so is conservative. Setting svdcut=None or svdcut=0 leaves the covariance matrix unchanged. Default is 1e-15.

Methods:
__call__(f, nopdf=False, **kargs)

Estimate expectation value of function f(p).

Uses module vegas to estimate the integral of f(p) multiplied by the probability density function associated with g (i.e., pdf(p)). The probability density function is omitted if nopdf=True (default setting is False).

Parameters
  • f (function) – Function f(p) to integrate. Integral is the expectation value of the function with respect to the distribution. The function can return a number, an array of numbers, or a dictionary whose values are numbers or arrays of numbers.

  • nopdf (bool) – If True drop the probability density function from the integrand (so no longer an expectation value). This is useful if one wants to use the optimized integrator for something other than a standard expectation value (e.g., an expectation value with a non-Gaussian PDF). Default is False.

All other keyword arguments are passed on to a vegas integrator; see the vegas documentation for further information.

Other Objects and Functions

class vegas.RAvg

Running average of scalar-valued Monte Carlo estimates.

This class accumulates independent Monte Carlo estimates (e.g., of an integral) and combines them into a single average. It is derived from gvar.GVar (from the gvar module if it is present) and represents a Gaussian random variable.

Different estimates are weighted by their inverse variances if parameter weight=True; otherwise straight, unweighted averages are used.

Attributes and methods:
mean

The mean value of the weighted average.

sdev

The standard deviation of the weighted average.

chi2

chi**2 of weighted average.

dof

Number of degrees of freedom in weighted average.

Q

Q or p-value of weighted average’s chi**2.

itn_results

A list of the results from each iteration.

add(g)

Add estimate g to the running average.

summary(weighted=None)

Assemble summary of results, iteration-by-iteration, into a string.

Parameters

weighted (bool) – Display weighted averages of results from different iterations if True; otherwise show unweighted averages. Default behavior is determined by vegas.

class vegas.RAvgArray

Running average of array-valued Monte Carlo estimates.

This class accumulates independent arrays of Monte Carlo estimates (e.g., of an integral) and combines them into an array of averages. It is derived from numpy.ndarray. The array elements are gvar.GVars (from the gvar module if present) and represent Gaussian random variables.

Different estimates are weighted by their inverse covariance matrices if parameter weight=True; otherwise straight, unweighted averages are used.

Attributes and methods:
chi2

chi**2 of weighted average.

dof

Number of degrees of freedom in weighted average.

Q

Q or p-value of weighted average’s chi**2.

itn_results

A list of the results from each iteration.

add(g)

Add estimate g to the running average.

summary(extended=False, weighted=None)

Assemble summary of results, iteration-by-iteration, into a string.

Parameters
  • extended (bool) – Include a table of final averages for every component of the integrand if True. Default is False.

  • weighted (bool) – Display weighted averages of results from different iterations if True; otherwise show unweighted averages. Default behavior is determined by vegas.

class vegas.RAvgDict

Running average of dictionary-valued Monte Carlo estimates.

This class accumulates independent dictionaries of Monte Carlo estimates (e.g., of an integral) and combines them into a dictionary of averages. It is derived from gvar.BufferDict. The dictionary values are gvar.GVars or arrays of gvar.GVars.

Different estimates are weighted by their inverse covariance matrices if parameter weight=True; otherwise straight, unweighted averages are used.

Attributes and methods:
chi2

chi**2 of weighted average.

dof

Number of degrees of freedom in weighted average.

Q

Q or p-value of weighted average’s chi**2.

itn_results

A list of the results from each iteration.

add(g)
summary(extended=False, weighted=None)

Assemble summary of results, iteration-by-iteration, into a string.

Parameters
  • extended (bool) – Include a table of final averages for every component of the integrand if True. Default is False.

  • weighted (bool) – Display weighted averages of results from different iterations if True; otherwise show unweighted averages. Default behavior is determined by vegas.

vegas.batchintegrand()

Decorator for batch integrand functions.

Applying vegas.batchintegrand() to a function fcn repackages the function in a format that vegas can understand. Appropriate functions take a numpy array of integration points x[i, d] as an argument, where i=0... labels the integration point and d=0... labels direction, and return an array f[i] of integrand values (or arrays of integrand values) for the corresponding points. The meaning of fcn(x) is unchanged by the decorator.

An example is

import vegas
import numpy as np

@vegas.batchintegrand
def f(x):
    return np.exp(-x[:, 0] - x[:, 1])

for the two-dimensional integrand \exp(-x_0 - x_1).

This decorator provides an alternative to deriving an integrand class from vegas.BatchIntegrand.

class vegas.BatchIntegrand

Base class for classes providing batch integrands.

A class derived from vegas.BatchIntegrand will normally provide a __call__(self, x) method that returns an array f where:

x[i, d] is a contiguous numpy array where i=0... labels different integrtion points and d=0... labels different directions in the integration space.

f[i] is a contiguous array containing the integrand values corresponding to the integration points x[i, :]. f[i] is either a number, for a single integrand, or an array (of any shape) for multiple integrands (i.e., an array-valued integrand).

An example is

import vegas
import numpy as np

class batchf(vegas.BatchIntegrand):
    def __call__(x):
        return np.exp(-x[:, 0] - x[:, 1])

f = batchf()      # the integrand

for the two-dimensional integrand \exp(-x_0 - x_1).

Deriving from vegas.BatchIntegrand is the easiest way to construct integrands in Cython, and gives the fastest results.