IntU.jl is a Julia package for the symbolic calculation of integrals over
the Haar measure of classical compact groups (
For detailed documentation, please visit iitis.github.io/IntU.jl.
To introduce the main functionality of IntU, consider the problem of averaging
IntU provides an exact analytic result instantly, even for symbolic dimensions, using a simple unified interface: integrate(expr, measure). It supports matrix-valued expressions and provides the @integrate macro for intuitive symbolic integration.
The @integrate macro implicitly identifies random matrices based on the measure:
-
dU,dSU$\rightarrow$ U(Unitary) -
dO$\rightarrow$ O(Orthogonal) -
dSp$\rightarrow$ Sp(Symplectic) -
dPerm$\rightarrow$ P(Permutation) -
dCPerm$\rightarrow$ Y(Centered Permutation) -
dCOE,dCSE$\rightarrow$ S(Circular Orthogonal/Symplectic) -
dPsi$\rightarrow$ psi(Pure State) -
dDiagUnitary$\rightarrow$ D(Diagonal Unitary) -
dStiefel$\rightarrow$ V(Stiefel Manifold) - Ginibre Ensembles
$\rightarrow$ G
Unknown symbols (like A, B, d) are automatically treated as constants or dimensions.
Note
Symbol Scope and Redefinition: The @integrate macro manages a persistent symbolic state in your session. If you use a symbol (e.g., U) as a random matrix in one call and then use it in a different context (e.g., as a constant in an orthogonal integral), the macro automatically re-declares and re-binds the symbol to the correct type and measure context. This "Safety Rebind" mechanism prevents silent math errors when running examples in sequence.
using IntU, Symbolics, LinearAlgebra
# Integrate the matrix expression U * U'
# The @integrate macro knows 'U' is the random matrix for measure dU(2)
res = @integrate U * U' dU(2)
# Output: [1.0 0.0; 0.0 1.0] (2x2 Identity Matrix)using IntU
# Compute the integral of |U_{1,1}|^2
@integrate abs(U[1, 1])^2 dU(d)
# Output: 1 / dFor more complex moments, such as
@integrate abs(U[1, 1])^2 * abs(U[1, 2])^2 dU(d)
# Output: 1 / (d * (1 + d))IntU implements Weingarten calculus for the Unitary, Orthogonal, and Symplectic groups, as well as Haar-random pure states.
Unitary matrices dU and integrate.
# 4-th moment of a diagonal entry
@integrate abs(U[1, 1])^4 dU(d)
# Output: 2 / (d * (1 + d))The Special Unitary group dSU.
@integrate abs(U[1, 1])^2 dSU(d)
# Output: 1/dOrthogonal matrices dO measure.
@integrate O[1, 1]^4 dO(d)
# Output: 3 / (d * (2 + d))Symplectic matrices dSp.
@integrate abs(Sp[1, 1])^2 dSp(d)
# Output: 1 / dGinibre ensembles consist of non-Hermitian matrices with i.i.d. Gaussian entries.
- GinUE (Complex Ginibre Ensemble): i.i.d. complex Gaussian entries. Use
dGinUE. - GinOE (Real Ginibre Ensemble): i.i.d. real Gaussian entries. Use
dGinOE. - GinSE (Symplectic Ginibre Ensemble): i.i.d. quaternionic Gaussian entries. Use
dGinSE.
# E[Tr(G G')] = d^2
@integrate tr(G * G') dGinUE(d)IntU also supports Circular Ensembles (CUE, COE, CSE) which are commonly used in random matrix theory.
-
CUE (Circular Unitary Ensemble): Equivalent to the Haar measure on
$U(d)$ . UsedCUE. -
COE (Circular Orthogonal Ensemble): Ensemble of symmetric unitary matrices. Use
dCOE. -
CSE (Circular Symplectic Ensemble): Ensemble of self-dual unitary matrices of even dimension
$2n$ . UsedCSE.
# COE moment E[|S_{1,1}|^2]
@integrate abs(S[1, 1])^2 dCOE(d)
# Output: 2 / (d + 1)IntU can integrate polynomial functions of the components of a Haar-random pure
state vector
# Average of |ψ_1|^2
@integrate abs(psi[1, 1])^2 dPsi(d)
# Output: 1 / dIntU supports integration over the Symmetric group
-
dPerm: Integration over the set of
$d \times d$ permutation matrices.
# E[P_11]
@integrate P[1, 1] dPerm(d)
# Output: 1 / dIntU also supports centered permutation matrices
- dCPerm: Integration over centered permutation matrices.
# E[Y_11^2]
@integrate Y[1, 1]^2 dCPerm(d)
# Output: (d - 1) / d^2IntU supports integration over the group of diagonal unitary matrices, which corresponds to independent phase averaging for each diagonal entry.
# E[|D_11|^2]
@integrate abs(D[1, 1])^2 dDiagUnitary(d)
# Output: 1IntU supports index-free notation for integrating traces of products of random matrices, which is often more convenient for quantum information tasks.
# Compute ∫ tr(U A U† B) dU
@integrate tr(U * A * U' * B) dU(d)
# Output: (tr(A)*tr(B)) / dIntU supports calculating HCIZ integrals of the form
using IntU, LinearAlgebra
A = diagm([1.0, 2.0])
B = diagm([0.5, 1.5])
# Compute ∫ dU exp(Tr(A U B U'))
hciz(A, B)
# Output: 20.9329...IntU supports integration over the Stiefel manifold
# E[|V_{1,1}|^2]
@integrate abs(V[1, 1])^2 dStiefel(d, 2)
# Output: 1 / dIntU.jl provides a bridge to ITensors.jl for symbolic integration of tensor networks.
using IntU, ITensors
i, j = Index(2), Index(2)
U_it = randomITensor(i, j)
U = ITensorUnitary(U_it; out_indices=[i], in_indices=[j])
A = randomITensor(j, i)
res = integrate([U, A], dU(2))IntU.jl comes with a comprehensive set of examples and benchmarks.
You can run all examples at once using the provided script:
sh examples/runexamples.shTo run a single example, pass the example number to the script. For instance, to run 01_basics.jl:
sh examples/runexamples.sh 01To run all benchmarks:
sh benchmarks/runbenchmarks.shTo run a single benchmark, pass the number to the script:
sh benchmarks/runbenchmarks.sh 01IntU is tested with Julia 1.11 or later. Installation can be done through the Pkg REPL:
import Pkg; Pkg.add(url="https://github.com/iitis/IntU.jl")If you use IntU.jl in your research, please cite:
@misc{intu2024,
author = {Pawela, Łukasz and Puchała, Zbigniew},
title = {IntU.jl: Symbolic integration over the Haar measure of classical compact groups},
year = {2024},
publisher = {GitHub},
journal = {GitHub repository},
howpublished = {\url{https://github.com/iitis/IntU.jl}}
}