From 949421f4eedbdf0493a6080a2ecd4f862f051c68 Mon Sep 17 00:00:00 2001 From: "Steven G. Johnson" Date: Fri, 17 May 2013 13:07:19 -0400 Subject: [PATCH] added quadgk integration routine --- base/exports.jl | 4 + base/quadgk.jl | 344 ++++++++++++++++++++++++++++++++++++++++++++ base/sysimg.jl | 4 + doc/helpdb.jl | 21 +++ doc/stdlib/base.rst | 54 +++++++ test/math.jl | 8 ++ 6 files changed, 435 insertions(+) create mode 100644 base/quadgk.jl diff --git a/base/exports.jl b/base/exports.jl index e98b677a94b17..0de7d3f71e37d 100644 --- a/base/exports.jl +++ b/base/exports.jl @@ -11,6 +11,7 @@ export Math, MPFR, GMP, + QuadGK, Sort, Test, Pkg, @@ -862,6 +863,9 @@ export plan_dct!, plan_idct!, +# numerical integration + quadgk, + # iteration start, done, diff --git a/base/quadgk.jl b/base/quadgk.jl new file mode 100644 index 0000000000000..685ac45d80ac6 --- /dev/null +++ b/base/quadgk.jl @@ -0,0 +1,344 @@ +module QuadGK +export gauss, kronrod, quadgk +using Base.Collections +import Base.isless, Base.Sort.Reverse + +# Adaptive Gauss-Kronrod quadrature routines (arbitrary precision), +# written and contributed to Julia by Steven G. Johnson, 2013. + +########################################################################### + +# cache of (T,n) -> (x,w,gw) Kronrod rules, to avoid recomputing them +# unnecessarily for repeated integration. We initialize it with the +# default n=7 rule for double-precision calculations. +const rulecache = (Any=>Any)[ (Float64,7) => # precomputed in 100-bit arith. + ([-9.9145537112081263920685469752598e-01, + -9.4910791234275852452618968404809e-01, + -8.6486442335976907278971278864098e-01, + -7.415311855993944398638647732811e-01, + -5.8608723546769113029414483825842e-01, + -4.0584515137739716690660641207707e-01, + -2.0778495500789846760068940377309e-01, + 0.0], + [2.2935322010529224963732008059913e-02, + 6.3092092629978553290700663189093e-02, + 1.0479001032225018383987632254189e-01, + 1.4065325971552591874518959051021e-01, + 1.6900472663926790282658342659795e-01, + 1.9035057806478540991325640242055e-01, + 2.0443294007529889241416199923466e-01, + 2.0948214108472782801299917489173e-01], + [1.2948496616886969327061143267787e-01, + 2.797053914892766679014677714229e-01, + 3.8183005050511894495036977548818e-01, + 4.1795918367346938775510204081658e-01]) ] + +# integration segment (a,b), estimated integral I, and estimated error E +immutable Segment + a::Number + b::Number + I::Number + E::Number +end +isless(i::Segment, j::Segment) = isless(i.E, j.E) + +# Internal routine: integrate f(x) over a sequence of line segments +# evaluate the integration rule and return a Segment. +function evalrule(f, a,b, x,w,gw) + # Ik and Ig are integrals via Kronrod and Gauss rules, respectively + Ik = Ig = zero(promote_type(typeof(a),eltype(w))) + s = convert(eltype(x), 0.5) * (b-a) + n1 = 1 - (length(x) & 1) # 0 if even order, 1 if odd order + for i = 1:length(gw)-n1 + fg = f(a + (1+x[2i])*s) + f(a + (1-x[2i])*s) + fk = f(a + (1+x[2i-1])*s) + f(a + (1-x[2i-1])*s) + Ig += fg * gw[i] + Ik += fg * w[2i] + fk * w[2i-1] + end + if n1 == 0 # even: Gauss rule does not include x == 0 + Ik += f(a + s) * w[end] + else # odd: don't count x==0 twice in Gauss rule + f0 = f(a + s) + Ig += f0 * gw[end] + Ik += f0 * w[end] + + (f(a + (1+x[end-1])*s) + f(a + (1-x[end-1])*s)) * w[end-1] + end + Ik *= s + Ig *= s + E = abs(Ik - Ig) + if isnan(E) || isinf(E) + throw(DomainError()) + end + return Segment(a, b, Ik, E) +end + +rulekey(::Type{BigFloat}, n) = (BigFloat, get_bigfloat_precision(), n) +rulekey(T,n) = (T,n) + +# Internal routine: integrate f over the union of the open intervals +# (s[1],s[2]), (s[2],s[3]), ..., (s[end-1],s[end]), using h-adaptive +# integration with the order-n Kronrod rule and weights of type Tw, +# with absolute tolerance abstol and relative tolerance reltol, +# with maxevals an approximate maximum number of f evaluations. +function do_quadgk{Tw}(f, s, n, ::Type{Tw}, abstol, reltol, maxevals) + + if eltype(s) <: Real # check for infinite or semi-infinite intervals + s1 = s[1]; s2 = s[end]; inf1 = isinf(s1); inf2 = isinf(s2) + if inf1 || inf2 + if inf1 && inf2 # x = t/(1-t^2) coordinate transformation + return do_quadgk(t -> begin t2 = t*t; den = 1 / (1 - t2); + f(t*den) * (1+t2)*den*den; end, + map(x -> isinf(x) ? copysign(one(x), x) : 2x / (1+hypot(1,2x)), s), + n, Tw, abstol, reltol, maxevals) + end + s0,si = inf1 ? (s2,s1) : (s1,s2) + if si < 0 # x = s0 - t/(1-t) + return do_quadgk(t -> begin den = 1 / (1 - t); + f(s0 - t*den) * den*den; end, + reverse!(map(x -> 1 / (1 + 1 / (s0 - x)), s)), + n, Tw, abstol, reltol, maxevals) + else # x = s0 + t/(1-t) + return do_quadgk(t -> begin den = 1 / (1 - t); + f(s0 + t*den) * den*den; end, + map(x -> 1 / (1 + 1 / (x - s0)), s), + n, Tw, abstol, reltol, maxevals) + end + end + end + + key = rulekey(Tw,n) + x,w,gw = haskey(rulecache, key) ? rulecache[key] : + (rulecache[key] = kronrod(Tw, n)) + segs = Segment[] + for i in 1:length(s) - 1 + heappush!(segs, evalrule(f, s[i],s[i+1], x,w,gw), Reverse) + end + numevals = (2n+1) * length(segs) + I = segs[1].I + E = segs[1].E + for i in 2:length(segs) + I += segs[i].I + E += segs[i].E + end + # Pop the biggest-error segment and subdivide (h-adaptation) + # until convergence is achieved or maxevals is exceeded. + while E > abstol && E > reltol * abs(I) && numevals < maxevals + s = heappop!(segs, Reverse) + mid = (s.a + s.b) * 0.5 + s1 = evalrule(f, s.a, mid, x,w,gw) + s2 = evalrule(f, mid, s.b, x,w,gw) + heappush!(segs, s1, Reverse) + heappush!(segs, s2, Reverse) + I = (I - s.I) + s1.I + s2.I + E = (E - s.E) + s1.E + s2.E + numevals += 4n+2 + end + # re-sum (paranoia about accumulated roundoff) + I = segs[1].I + E = segs[1].E + for i in 2:length(segs) + I += segs[i].I + E += segs[i].E + end + return (I, E) +end + +# Gauss-Kronrod quadrature of f from a to b to c... + +function quadgk{T<:FloatingPoint}(f, a::T,b::T,c::T...; + abstol=zero(T), reltol=eps(T)*100, + maxevals=10^7, order=7) + do_quadgk(f, [a, b, c...], order, T, abstol, reltol, maxevals) +end + +function quadgk{T<:FloatingPoint}(f, a::Complex{T}, + b::Complex{T},c::Complex{T}...; + abstol=zero(T), reltol=eps(T)*100, + maxevals=10^7, order=7) + do_quadgk(f, [a, b, c...], order, T, abstol, reltol, maxevals) +end + +# generic version: determine precision from a combination of +# all the integration-segment endpoints +function quadgk(f, a, b, c...; kws...) + T = promote_type(typeof(float(a)), typeof(b)) + for i in 1:length(c) + T = promote_type(T, typeof(c[i])) + end + cT = T[ c[i] for i in 1:length(c) ] + quadgk(f, convert(T, a), convert(T, b), cT...; kws...) +end + +########################################################################### +# Gauss-Kronrod integration-weight computation for arbitrary floating-point +# types and precision, implemented based on the description in: +# +# Dirk P. Laurie, "Calculation of Gauss-Kronrod quadrature rules," +# Mathematics of Computation, vol. 66, no. 219, pp. 1133-1145 (1997). +# +# for the Kronrod rule, and for the Gauss rule from the description in +# +# Lloyd N. Trefethen and David Bau, Numerical Linear Algebra (SIAM, 1997). +# +# Arbitrary-precision eigenvalue (eignewt & eigpoly) and eigenvector +# (eigvec1) routines are written by SGJ, independent of the above sources. +# +# Since we only implement Gauss-Kronrod rules for the unit weight function, +# the diagonals of the Jacobi matrices are zero and certain things simplify +# compared to the general case of an arbitrary weight function. + +# Given a symmetric tridiagonal matrix H with H[i,i] = 0 and +# H[i-1,i] = H[i,i-1] = b[i-1], compute p(z) = det(z I - H) and its +# derivative p'(z), returning (p,p'). +function eigpoly(b,z,m=length(b)+1) + d1 = z + d1deriv = d2 = one(z); + d2deriv = zero(z); + for i = 2:m + b2 = b[i-1]^2 + d = z * d1 - b2 * d2; + dderiv = d1 + z * d1deriv - b2 * d2deriv; + d2 = d1; + d1 = d; + d2deriv = d1deriv; + d1deriv = dderiv; + end + return (d1, d1deriv) +end + +# compute the n smallest eigenvalues of the symmetric tridiagonal matrix H +# (defined from b as in eigpoly) using a Newton iteration +# on det(H - lambda I). Unlike eig, handles BigFloat. +function eignewt(b,m,n) + # get initial guess from eig on Float64 matrix + H = zeros(m, m) + for i in 2:m + H[i-1,i] = H[i,i-1] = b[i-1] + end + # TODO: use LAPACK routine for eigenvalues of tridiagonal matrix + # rather than generic O(n^3) routine + lambda0 = sort(eigvals(H)) + + lambda = Array(eltype(b), n) + for i = 1:n + lambda[i] = lambda0[i] + for k = 1:1000 + (p,pderiv) = eigpoly(b,lambda[i],m) + lambda[i] = (lamold = lambda[i]) - p / pderiv + if abs(lambda[i] - lamold) < 10 * eps(lambda[i]) * abs(lambda[i]) + break + end + end + # do one final Newton iteration for luck and profit: + (p,pderiv) = eigpoly(b,lambda[i],m) + lambda[i] = lambda[i] - p / pderiv + end + return lambda +end + +# given an eigenvalue z and the matrix H(b) from above, return +# the corresponding eigenvector, normalized to 1. +function eigvec1(b,z::Number,m=length(b)+1) + # "cheat" and use the fact that our eigenvector v must have a + # nonzero first entries (since it is a quadrature weight), so we + # can set v[1] = 1 to solve for the rest of the components:. + v = Array(eltype(b), m) + v[1] = 1 + if m > 1 + s = v[1] + v[2] = z * v[1] / b[1] + s += v[2]^2 + for i = 3:m + v[i] = - (b[i-2]*v[i-2] - z*v[i-1]) / b[i-1] + s += v[i]^2 + end + scale!(v, 1 / sqrt(s)) + end + return v +end + +# return (x, w) pair of N quadrature points x[i] and weights w[i] to +# integrate functions on the interval (-1, 1). i.e. dot(f(x), w) +# approximates the integral. Uses the method described in Trefethen & +# Bau, Numerical Linear Algebra, to find the N-point Gaussian quadrature +function gauss{T<:FloatingPoint}(::Type{T}, N::Integer) + if N < 1 + throw(ArgumentError("Gauss rules require positive order")) + end + o = one(T) + b = T[ n / sqrt(4n^2 - o) for n = 1:N-1 ] + x = eignewt(b,N,N) + w = T[ 2*eigvec1(b,x[i])[1]^2 for i = 1:N ] + return (x, w) +end + +# Compute 2n+1 Kronrod points x and weights w based on the description in +# Laurie (1997), appendix A, simplified for a=0, for integrating on [-1,1]. +# Since the rule is symmetric only returns the n+1 points with x <= 0. +# Also computes the embedded n-point Gauss quadrature weights gw (again +# for x <= 0), corresponding to the points x[2:2:end]. Returns (x,w,wg). +function kronrod{T<:FloatingPoint}(::Type{T}, n::Integer) + if n < 1 + throw(ArgumentError("Kronrod rules require positive order")) + end + o = one(T) + b = zeros(T, 2n+1) + b[1] = 2*o + for j = 1:div(3n+1,2) + b[j+1] = j^2 / (4j^2 - o) + end + s = zeros(T, div(n,2) + 2) + t = zeros(T, div(n,2) + 2) + t[2] = b[n+2] + for m = 0:n-2 + u = zero(T) + for k = div(m+1,2):-1:0 + l = m - k + 1 + k1 = k + n + 2 + u += b[k1]*s[k+1] - b[l]*s[k+2] + s[k+2] = u + end + s,t = t,s + end + for j = div(n,2):-1:0 + s[j+2] = s[j+1] + end + for m = n-1:2n-3 + u = zero(T) + for k = m+1-n:div(m-1,2) + l = m - k + 1 + j = n - l + k1 = k + n + 2 + u -= b[k1]*s[j+2] - b[l]*s[j+3] + s[j+2] = u + end + k = div(m+1,2) + if 2k != m + j = n - (m - k + 2) + b[k+n+2] = s[j+2] / s[j+3] + end + s,t = t,s + end + for j = 1:2n + b[j] = sqrt(b[j+1]) + end + + # get negative quadrature points x + x = eignewt(b,2n+1,n+1) # x <= 0 + + # get quadrature weights + w = T[ 2*eigvec1(b,x[i],2n+1)[1]^2 for i in 1:n+1 ] + + # Get embedded Gauss rule from even-indexed points, using + # the method described in Trefethen and Bau. + for j = 1:n-1 + b[j] = j / sqrt(4j^2 - o) + end + gw = T[ 2*eigvec1(b,x[i],n)[1]^2 for i = 2:2:n+1 ] + + return (x, w, gw) +end + +########################################################################### + +end # module QuadGK diff --git a/base/sysimg.jl b/base/sysimg.jl index 3d93508ea2d8c..d9b9e37eb668e 100644 --- a/base/sysimg.jl +++ b/base/sysimg.jl @@ -176,6 +176,10 @@ importall .GMP include("mpfr.jl") importall .MPFR +# Numerical integration +include("quadgk.jl") +importall .QuadGK + # deprecated functions include("deprecated.jl") diff --git a/doc/helpdb.jl b/doc/helpdb.jl index b9541ad79d7d5..388ba8a580e3c 100644 --- a/doc/helpdb.jl +++ b/doc/helpdb.jl @@ -4180,6 +4180,27 @@ "), +("Numerical Integration","Base","quadgk","quadgk(f, a,b,c...; options) + + Numerically integrate the function \"f(x)\" from \"a\" to \"b\", + and optionally over additional intervals \"b\" to \"c\" and so on. + Keyword options include a relative error tolerance \"reltol\" (defaults + to \"100*eps\"), an absolute error tolerance \"abstol\" (defaults + to 0), a maximum number of function evaluations \"maxevals\" (defaults + to \"10^7\"), and the \"order\" of the integration rule (defaults to 7). + + Returns a pair \"(I,E)\" of the estimated integral \"I\" and an + estimated upper bound on the absolute error \"E\". + + Complex-valued functions are supported, and the endpoints \"a\" etcetera + can also be complex (in which case the integral is performed over + straight-line segments in the complex plane). If the endpoints + are \"BigFloat\", then the integration will be performed in that + precision as well (note: it is advisable to increase the integration + \"order\" in rough proportion to the precision, for smooth integrands). + +"), + ("Parallel Computing","Base","addprocs","addprocs(n) Add processes on the local machine. Can be used to take advantage diff --git a/doc/stdlib/base.rst b/doc/stdlib/base.rst index 3bc99a36da5ea..8e0f61eb1e1f3 100644 --- a/doc/stdlib/base.rst +++ b/doc/stdlib/base.rst @@ -2777,6 +2777,60 @@ FFT functions in Julia are largely implemented by calling functions from `FFTW < Compute the cross-correlation of two vectors. +Numerical Integration +----------------- + +Although several external packages are available for numeric integration +and solution of ordinary differential equations, we also provide +some built-in integration support in Julia. + +.. function:: quadgk(f, a,b,c...; reltol=eps*100, abstol=0, maxevals=10^7, order=7) + + Numerically integrate the function ``f(x)`` from ``a`` to ``b``, + and optionally over additional intervals ``b`` to ``c`` and so on. + Keyword options include a relative error tolerance ``reltol`` (defaults + to ``100*eps`` in the precision of the endpoints), an absolute error + tolerance ``abstol`` (defaults to 0), a maximum number of function + evaluations ``maxevals`` (defaults to ``10^7``), and the ``order`` + of the integration rule (defaults to 7). + + Returns a pair ``(I,E)`` of the estimated integral ``I`` and an + estimated upper bound on the absolute error ``E``. If ``maxevals`` + is not exceeded then either ``E <= abstol`` or ``E <= + reltol*abs(I)`` will hold. (Note that it is useful to specify a + positive ``abstol`` in cases where ``abs(I)`` may be zero.) + + Complex-valued functions are supported, and the endpoints ``a`` etcetera + can also be complex (in which case the integral is performed over + straight-line segments in the complex plane). If the endpoints + are ``BigFloat``, then the integration will be performed in ``BigFloat`` + precision as well (note: it is advisable to increase the integration + ``order`` in rough proportion to the precision, for smooth integrands). + More generally, the precision is set by the precision of the integration + endpoints (promoted to floating-point types). + + The algorithm is an adaptive Gauss-Kronrod integration technique: + the integral in each interval is estimated using a Kronrod rule + (``2*order+1`` points) and the error is estimated using an embedded + Gauss rule (``order`` points). The interval with the largest + error is then subdivided into two intervals and the process is repeated + until the desired error tolerance is achieved. + + These quadrature rules work best for smooth functions within each + interval, so if your function has a known discontinuity or other + singularity, it is best to subdivide your interval to put the + singularity at an endpoint. For example, if ``f`` has a discontinuity + at ``x=0.7`` and you want to integrate from 0 to 1, you should use + ``quadgk(f, 0,0.7,1)`` to subdivide the interval at the point of + discontinuity. The integrand is never evaluated exactly at the endpoints + of the intervals, so it is possible to integrate functions that diverge + at the endpoints as long as the singularity is integrable (for example, + a ``log(x)`` or ``1/sqrt(x)`` singularity). + + For real-valued endpoints, the starting and/or ending points may be + infinite. (A coordinate transformation is performed internally to + map the infinite interval to a finite one.) + Parallel Computing ------------------ diff --git a/test/math.jl b/test/math.jl index f16902f35b32e..691771882441e 100644 --- a/test/math.jl +++ b/test/math.jl @@ -121,3 +121,11 @@ end @test_approx_eq zeta(0) -0.5 @test_approx_eq zeta(2) pi^2/6 @test_approx_eq zeta(4) pi^4/90 + +@test_approx_eq quadgk(cos, 0,0.7,1)[1] sin(1) +@test_approx_eq quadgk(x -> exp(im*x), 0,0.7,1)[1] (exp(1im)-1)/im +@test_approx_eq quadgk(x -> exp(im*x), 0,1im)[1] -1im*expm1(-1) +@test_approx_eq_eps quadgk(cos, 0,BigFloat(1),order=40)[1] sin(BigFloat(1)) 1000*eps(BigFloat) +@test_approx_eq quadgk(x -> exp(-x), 0,0.7,Inf)[1] 1.0 +@test_approx_eq quadgk(x -> exp(x), -Inf,0)[1] 1.0 +@test_approx_eq quadgk(x -> exp(-x^2), -Inf,Inf)[1] sqrt(pi)