Mitsuba Renderer  0.5.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
mitsuba::NDIntegrator Class Reference

Adaptively computes the integral of a multidimensional function using either a Gauss-Kronod (1D) or a Genz-Malik (>1D) cubature rule. More...

#include <mitsuba/core/quad.h>

Public Types

enum  EResult { ESuccess = 0, EFailure = 1 }
 
typedef boost::function< void(const
Float *, Float *)> 
Integrand
 
typedef boost::function< void(size_t,
const Float *, Float *)> 
VectorizedIntegrand
 

Public Member Functions

 NDIntegrator (size_t fDim, size_t dim, size_t maxEvals, Float absError=0, Float relError=0)
 
EResult integrate (const Integrand &f, const Float *min, const Float *max, Float *result, Float *error, size_t *evals=NULL) const
 Integrate the function f over the rectangular domain bounded by min and max. More...
 
EResult integrateVectorized (const VectorizedIntegrand &f, const Float *min, const Float *max, Float *result, Float *error, size_t *evals=NULL) const
 Integrate the function f over the rectangular domain bounded by min and max. More...
 

Protected Attributes

size_t m_fdim
 
size_t m_dim
 
size_t m_maxEvals
 
Float m_absError
 
Float m_relError
 

Detailed Description

Adaptively computes the integral of a multidimensional function using either a Gauss-Kronod (1D) or a Genz-Malik (>1D) cubature rule.

This class is a C++ wrapper around the cubature code by Steven G. Johnson (http://ab-initio.mit.edu/wiki/index.php/Cubature)

The original implementation is based on algorithms proposed in

A. C. Genz and A. A. Malik, "An adaptive algorithm for numeric integration over an N-dimensional rectangular region," J. Comput. Appl. Math. 6 (4), 295–302 (1980).

and

J. Berntsen, T. O. Espelid, and A. Genz, "An adaptive algorithm for the approximate calculation of multiple integrals," ACM Trans. Math. Soft. 17 (4), 437–451 (1991).

Member Typedef Documentation

typedef boost::function<void (const Float *, Float *)> mitsuba::NDIntegrator::Integrand
typedef boost::function<void (size_t, const Float *, Float *)> mitsuba::NDIntegrator::VectorizedIntegrand

Member Enumeration Documentation

Enumerator
ESuccess 
EFailure 

Constructor & Destructor Documentation

mitsuba::NDIntegrator::NDIntegrator ( size_t  fDim,
size_t  dim,
size_t  maxEvals,
Float  absError = 0,
Float  relError = 0 
)

Initialize the Cubature integration scheme

Parameters
fDimNumber of integrands (i.e. dimensions of the image space)
nDimNumber of integration dimensions (i.e. dimensions of the function domain)
maxEvalsMaximum number of function evaluations (0 means no limit). The error bounds will likely be exceeded when the integration is forced to stop prematurely. Note: the actual number of evaluations may somewhat exceed this value.
absErrorAbsolute error requirement (0 to disable)
relErrorRelative error requirement (0 to disable)

Member Function Documentation

EResult mitsuba::NDIntegrator::integrate ( const Integrand f,
const Float min,
const Float max,
Float result,
Float error,
size_t *  evals = NULL 
) const

Integrate the function f over the rectangular domain bounded by min and max.

The supplied function should have the interface

void integrand(const Float *in, Float *out);

The input array in consists of one set of input parameters having dim entries. The function is expected to store the results of the evaluation into the out array using fDim entries.

EResult mitsuba::NDIntegrator::integrateVectorized ( const VectorizedIntegrand f,
const Float min,
const Float max,
Float result,
Float error,
size_t *  evals = NULL 
) const

Integrate the function f over the rectangular domain bounded by min and max.

This function implements a vectorized version of the above integration function, which is more efficient by evaluating the integrant in `batches'. The supplied function should have the interface

void integrand(int numPoints, const Float *in, Float *out);

Note that in in is not a single point, but an array of numPoints points (length numPoints x dim), and upon return the values of all fDim integrands at all numPoints points should be stored in out (length fDim x numPoints). In particular, out[i*dim + j] is the j-th coordinate of the i-th point, and the k-th function evaluation (k<fDim) for the i-th point is returned in out[k*npt + i]. The size of numPoints will vary with the dimensionality of the problem; higher-dimensional problems will have (exponentially) larger numbers, allowing for the possibility of more parallelism. Currently, numPoints starts at 15 in 1d, 17 in 2d, and 33 in 3d, but as the integrator calls your integrand more and more times the value will grow. e.g. if you end up requiring several thousand points in total, numPoints may grow to several hundred.

Member Data Documentation

Float mitsuba::NDIntegrator::m_absError
protected
size_t mitsuba::NDIntegrator::m_dim
protected
size_t mitsuba::NDIntegrator::m_fdim
protected
size_t mitsuba::NDIntegrator::m_maxEvals
protected
Float mitsuba::NDIntegrator::m_relError
protected

The documentation for this class was generated from the following file: