Integral Solver Algorithms
The following algorithms are available:
QuadGKJL
: Uses QuadGK.jl, which supports one-dimensional integration of scalar and array-valued integrands with in-place or batched forms. Integrands that are both in-place and batched are implemented in the wrapper but are not supported under the hood.HCubatureJL
: Uses HCubature.jl, which supports scalar and array-valued integrands and works best in low dimensions, e.g. ≤ 8. In-place integrands are implemented in the wrapper but are not supported under the hood. Batching is not supported.VEGAS
: Uses MonteCarloIntegration.jl, which requires scalar,Float64
-valued integrands and works in any number of dimensions.VEGASMC
: Uses MCIntegration.jl. Requiresusing MCIntegration
. Doesn't support batching.CubatureJLh
: h-Cubature from Cubature.jl. Requiresusing Cubature
.CubatureJLp
: p-Cubature from Cubature.jl. Requiresusing Cubature
.CubaVegas
: Vegas from Cuba.jl. Requiresusing Cuba
.CubaSUAVE
: SUAVE from Cuba.jl. Requiresusing Cuba
.CubaDivonne
: Divonne from Cuba.jl. Requiresusing Cuba
. Works only for>1
-dimensional integrations.CubaCuhre
: Cuhre from Cuba.jl. Requiresusing Cuba
. Works only for>1
-dimensional integrations.GaussLegendre
: Performs fixed-order Gauss-Legendre quadrature. Requiresusing FastGaussQuadrature
.QuadratureRule
: Accepts a user-defined function that returns nodes and weights.ArblibJL
: real- and complex-valued univariate integration of holomorphic and meromorphic functions from Arblib.jl. Requiresusing Arblib
.
Integrals.QuadGKJL
— TypeQuadGKJL(; order = 7, norm = norm, buffer = nothing)
One-dimensional Gauss-Kronrod integration from QuadGK.jl. This method also takes the optional arguments order
and norm
. Which are the order of the integration rule and the norm for calculating the error, respectively. Lastly, the buffer
keyword, if set (e.g. buffer=true
), will allocate a buffer to reuse for multiple integrals and may require evaluating the integrand unless an integrand_prototype
is provided. Unlike the segbuf
keyword to quadgk
, you do not allocate the buffer as this is handled automatically.
References
@article{laurie1997calculation,
title={Calculation of Gauss-Kronrod quadrature rules},
author={Laurie, Dirk},
journal={Mathematics of Computation},
volume={66},
number={219},
pages={1133--1145},
year={1997}
}
Integrals.HCubatureJL
— TypeHCubatureJL(; norm=norm, initdiv=1, buffer = nothing)
Multidimensional "h-adaptive" integration from HCubature.jl. This method also takes the optional arguments initdiv
and norm
. Which are the initial number of segments each dimension of the integration domain is divided into, and the norm for calculating the error, respectively. Lastly, the buffer
keyword, if set (e.g. buffer=true
), will allocate a buffer to reuse for multiple integrals and may require evaluating the integrand unless an integrand_prototype
is provided. Unlike the buffer
keyword to hcubature/hquadrature
, you do not allocate the buffer as this is handled automatically.
References
@article{genz1980remarks,
title={Remarks on algorithm 006: An adaptive algorithm for numerical integration over an N-dimensional rectangular region},
author={Genz, Alan C and Malik, Aftab Ahmad},
journal={Journal of Computational and Applied mathematics},
volume={6},
number={4},
pages={295--302},
year={1980},
publisher={Elsevier}
}
Integrals.CubatureJLp
— TypeCubatureJLp(; error_norm=Cubature.INDIVIDUAL)
Multidimensional p-adaptive integration from Cubature.jl. This method is based on repeatedly doubling the degree of the cubature rules, until convergence is achieved. The used cubature rule is a tensor product of Clenshaw–Curtis quadrature rules. error_norm
specifies the convergence criterion for vector valued integrands. Defaults to Cubature.INDIVIDUAL
, other options are Cubature.PAIRED
, Cubature.L1
, Cubature.L2
, or Cubature.LINF
.
Integrals.CubatureJLh
— TypeCubatureJLh(; error_norm=Cubature.INDIVIDUAL)
Multidimensional h-adaptive integration from Cubature.jl. error_norm
specifies the convergence criterion for vector valued integrands. Defaults to Cubature.INDIVIDUAL
, other options are Cubature.PAIRED
, Cubature.L1
, Cubature.L2
, or Cubature.LINF
.
References
@article{genz1980remarks,
title={Remarks on algorithm 006: An adaptive algorithm for numerical integration over an N-dimensional rectangular region},
author={Genz, Alan C and Malik, Aftab Ahmad},
journal={Journal of Computational and Applied mathematics},
volume={6},
number={4},
pages={295--302},
year={1980},
publisher={Elsevier}
}
Integrals.VEGAS
— TypeVEGAS(; nbins = 100, ncalls = 1000, debug=false, seed = nothing)
Multidimensional adaptive Monte Carlo integration from MonteCarloIntegration.jl. Importance sampling is used to reduce variance. This method also takes three optional arguments nbins
, ncalls
and debug
which are the initial number of bins each dimension of the integration domain is divided into, the number of function calls per iteration of the algorithm, and whether debug info should be printed, respectively.
Limitations
This algorithm can only integrate Float64
-valued functions
References
@article{lepage1978new,
title={A new algorithm for adaptive multidimensional integration},
author={Lepage, G Peter},
journal={Journal of Computational Physics},
volume={27},
number={2},
pages={192--203},
year={1978},
publisher={Elsevier}
}
Integrals.VEGASMC
— TypeVEGASMC(; kws...)
Markov-chain based Vegas algorithm from MCIntegration.jl
Refer to MCIntegration.integrate
for documentation on the keywords, which are passed directly to the solver with a set of defaults that works for conforming integrands.
Integrals.CubaVegas
— TypeCubaVegas()
Multidimensional adaptive Monte Carlo integration from Cuba.jl. Importance sampling is used to reduce variance.
References
@article{lepage1978new,
title={A new algorithm for adaptive multidimensional integration},
author={Lepage, G Peter},
journal={Journal of Computational Physics},
volume={27},
number={2},
pages={192--203},
year={1978},
publisher={Elsevier}
}
Integrals.CubaSUAVE
— TypeCubaSUAVE()
Multidimensional adaptive Monte Carlo integration from Cuba.jl. Suave stands for subregion-adaptive VEGAS. Importance sampling and subdivision are thus used to reduce variance.
References
@article{hahn2005cuba,
title={Cuba—a library for multidimensional numerical integration},
author={Hahn, Thomas},
journal={Computer Physics Communications},
volume={168},
number={2},
pages={78--95},
year={2005},
publisher={Elsevier}
}
Integrals.CubaDivonne
— TypeCubaDivonne()
Multidimensional adaptive Monte Carlo integration from Cuba.jl. Stratified sampling is used to reduce variance.
References
@article{friedman1981nested,
title={A nested partitioning procedure for numerical multiple integration},
author={Friedman, Jerome H and Wright, Margaret H},
journal={ACM Transactions on Mathematical Software (TOMS)},
volume={7},
number={1},
pages={76--92},
year={1981},
publisher={ACM New York, NY, USA}
}
Integrals.CubaCuhre
— TypeCubaCuhre()
Multidimensional h-adaptive integration from Cuba.jl.
References
@article{berntsen1991adaptive,
title={An adaptive algorithm for the approximate calculation of multiple integrals},
author={Berntsen, Jarle and Espelid, Terje O and Genz, Alan},
journal={ACM Transactions on Mathematical Software (TOMS)},
volume={17},
number={4},
pages={437--451},
year={1991},
publisher={ACM New York, NY, USA}
}
Integrals.GaussLegendre
— TypeGaussLegendre{C, N, W}
Struct for evaluating an integral via (composite) Gauss-Legendre quadrature. The field C
will be true
if subintervals > 1
, and false
otherwise.
The fields nodes::N
and weights::W
are defined by nodes, weights = gausslegendre(n)
for a given number of nodes n
.
The field subintervals::Int64 = 1
(with default value 1
) defines the number of intervals to partition the original interval of integration [a, b]
into, splitting it into [xⱼ, xⱼ₊₁]
for j = 1,…,subintervals
, where xⱼ = a + (j-1)h
and h = (b-a)/subintervals
. Gauss-Legendre quadrature is then applied on each subinterval. For example, if [a, b] = [-1, 1]
and subintervals = 2
, then Gauss-Legendre quadrature will be applied separately on [-1, 0]
and [0, 1]
, summing the two results.
Integrals.QuadratureRule
— TypeQuadratureRule(q; n=250)
Algorithm to construct and evaluate a quadrature rule q
of n
points computed from the inputs as x, w = q(n)
. It assumes the nodes and weights are for the standard interval [-1, 1]^d
in d
dimensions, and rescales the nodes to the specific hypercube being solved. The nodes x
may be scalars in 1d or vectors in arbitrary dimensions, and the weights w
must be scalar. The algorithm computes the quadrature rule sum(w .* f.(x))
and the caller must check that the result is converged with respect to n
.
Integrals.ArblibJL
— TypeArblibJL(; check_analytic=false, take_prec=false, warn_on_no_convergence=false, opts=C_NULL)
One-dimensional adaptive Gauss-Legendre integration using rigorous error bounds and precision ball arithmetic. Generally this assumes the integrand is holomorphic or meromorphic, which is the user's responsibility to verify. The result of the integral is not guaranteed to satisfy the requested tolerances, however the result is guaranteed to be within the error estimate.
Arblib.jl only supports integration of univariate real- and complex-valued functions with both inplace and out-of-place forms. See their documentation for additional details the algorithm arguments and on implementing high-precision integrands. Additionally, the error estimate is included in the return value of the integral, representing a ball.