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.
These are described in detail below.
The central component of the vegas package is the integrator class:
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 can be array-valued, in which case f(x) returns an array of values corresponding to different integrands. Also 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., an hour) because it allows you to monitor progress as it is being made (or not).
Parameters: |
|
---|
vegas.Integrator objects have attributes for each of these parameters. In addition they have the following methods:
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. Array-valued integrands are useful for 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).
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)
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. |
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.
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() 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. The number in a batch is controlled by parameter nhcube_batch.
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:
where and
are the limits of integration.
The grid specifies the transformation function at the points
for
:
Linear interpolation is used between those points. The Jacobian for this transformation is:
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.
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: |
|
---|
Number of dimensions.
Number of increments along each grid axis.
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.
The increment widths of the grid:
self.inc[d, i] = self.grid[d, i + 1] - self.grid[d, i]
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: |
|
---|
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: |
|
---|
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).
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).
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 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, float)
jac = numpy.empty(y.shape[0], float)
Parameters: |
|
---|
Display plots showing the current grid.
Parameters: |
|
---|
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).
Running average of 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.
The mean value of the weighted average.
The standard deviation of the weighted average.
chi**2 of weighted average.
Number of degrees of freedom in weighted average.
Q or p-value of weighted average’s chi**2.
A list of the results from each iteration.
Add estimate g to the running average.
Assemble summary of independent results into a string.
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.
chi**2 of weighted average.
Number of degrees of freedom in weighted average.
Q or p-value of weighted average’s chi**2.
A list of the results from each iteration.
Add estimate g to the running average.
Assemble summary of independent results into a string.
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, but the type of fcn is changed.
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 .
This decorator provides an alternative to deriving an integrand class from 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 .
Deriving from vegas.BatchIntegrand is the easiest way to construct integrands in Cython, and gives the fastest results.
Convert (batch) integrand into an MPI multiprocessor integrand.
Applying decorator vegas.MPIintegrand to a function repackages the function as a batch vegas integrand that can execute in parallel on multiple processors. 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 f[i,...] of integrand values) for the corresponding points.
An example is
import vegas
import numpy as np
@vegas.MPIintegrand
def f(x):
return np.exp(-x[:, 0] - x[:, 1])
for the two-dimensional integrand . Of
course, one could write f = vegas.MPIintegrand(f) instead of
using the decorator.
Message passing between processors uses MPI via Python module mpi4py, which must be installed in Python. To run an MPI integration code mpi-integral.py on 4 processors, for example, one might execute:
mpirun -np 4 python mpi-integral.py
Executing python mpi-integral.py, without the mpirun, causes it to run on a single processor, in more or less the same way an integral with a batch integrand runs.
An object of type vegas.MPIintegrand contains information about the MPI processes in the following attributes:
MPI intracommunicator — mpi4py.MPI.Intracomm object mpi4py.MPI.COMM_WORLD.
Number of processors used.
MPI rank of current process. Each process has a unique rank, ranging from 0 to nproc-1. The rank is used to make different processes do different things (for example, one generally wants only one of the processes to report out final results).
The random number see used to reset numpy.random.random in all the processes.
The implementation used here has the entire integration code run on every processor. It is only when evaluating the integrand that the processes do different things. This is efficient provided most of the time is spent evaluating the integrand, which, in any case, is the only situation where it might make sense to use multiple processors.
Note that vegas.MPIintegrand assumes that vegas.Integrator is using the default random number generator (numpy.random.random). If this is not the case, it is important to seed the other random number generator so that all processes use the same random numbers.
The approach used here to make vegas parallel is based on a strategy developed by R. Horgan and Q. Mason for the original Fortran version of vegas.