API Reference
Core Integration
IntegrateUnitary.integrate — Function
integrate(expr::SymbolicMatrix, measure)Integrate a SymbolicMatrix as a whole. Returns a matrix of results if the matrix dimensions are concrete integers. Supports rectangular matrices (e.g., pure states with size (d, 1) or Stiefel matrices with size (d, k)).
integrate(expr::SymbolicMatrixProduct, measure)Integrate a product of SymbolicMatrices. Returns a matrix of results if the dimension in measure is a concrete integer. It skips expansion and returns the product itself if it does not contain the integration variable.
integrate(expr::AbstractArray, measure)Performs element-wise integration of a matrix or array of expressions. Returns an array of the same shape containing integrated values.
integrate(expr, measure)Top-level integration function. It first checks the Pre-computed Integral Library for instant results. If not found, it calls fallback_integrate for the specific measure.
IntegrateUnitary.@integrate — Macro
@integrate expr measureSymbolically integrate an expression over a measure, automatically declaring variables. Heuristics are used to identify which symbols represent random matrices and dimensions.
IntegrateUnitary.evaluate — Function
evaluate(expr, dict)
evaluate(expr, pair)Shorthand for Symbolics.substitute. Useful for substituting symbolic dimensions with numeric values in integration results. Also handles substituting symbolic traces.
This function automatically handles removable singularities in fractions (e.g., $0/0$ forms). If a denominator evaluates to zero after substitution, the expression is first simplified to attempt resolving the singularity before completing the evaluation.
IntegrateUnitary.asymptotic — Function
asymptotic(expr, measure::AbstractMeasure, order=1)Returns the series expansion of the integral in powers of 1/d. Works for any measure that implements _reconstruct_symbolic.
asymptotic(ex, d, order=1)Generic asymptotic expansion of a rational function ex in powers of 1/d.
IntegrateUnitary.hciz — Function
hciz(A, B)
hciz(a::AbstractVector, b::AbstractVector)Computes the Harish-Chandra-Itzykson-Zuber (HCIZ) integral:
\[\int_{U(d)} dU e^{\text{Tr}(A U B U^\dagger)} = \left( \prod_{p=1}^{d-1} p! \right) \frac{\det(e^{a_i b_j})_{i,j=1}^d}{\Delta(a) \Delta(b)}\]
where a and b are eigenvalues of A and B, and \Delta is the Vandermonde determinant.
If A and B are matrices, their eigenvalues are extracted. Supporting:
- Numeric matrices (via
eigen) Matrix{Num}(symbolic diagonal or 2x2)SymbolicMatrixwith concrete integer dimension (generates symbolic eigenvaluesa_1, ..., a_d)
Note: This formula is sensitive to degenerate eigenvalues where the denominators become zero. In such cases, the limit should be taken. This implementation currently uses a small perturbation for numerical stability if exact degeneracies are detected in numeric input.
IntegrateUnitary.vandermonde_det — Function
vandermonde_det(v)Computes the Vandermonde determinant of a vector v:
\[\Delta(v) = \prod_{1 \le i < j \le d} (v_i - v_j)\]
[!NOTE]
integrate(expr, measure)is the universal entry point for all calculations in IntegrateUnitary.jl. It handles matrix-valued expressions and library lookups, and supports symbolic dimensions for a broad class of entry-wise and trace-polynomial workflows. Some paths require concrete integer dimensions (including|tr(U)|^(2k)fork > 1,hcizonSymbolicMatrixinputs, and direct matrix-valued integration ofSymbolicMatrix/SymbolicMatrixProductexpressions).
Measures
IntegrateUnitary.AbstractMeasure — Type
AbstractMeasureAbstract base type for all integration measures (Haar, Gaussian, Circular, etc.). Provides generic dispatch for element-wise integration of arrays and ambiguity resolution for SymbolicMatrix and SymbolicMatrixProduct.
IntegrateUnitary.measure_info — Function
measure_info(measure)Returns (subs_dict, matcher, dim, measure_type) for a given measure. Subtypes participate in the unified integration flow by defining _measure_tag. Measures with custom logic can override measure_info directly.
Unitary Group
IntegrateUnitary.dU — Function
dU(dim)Defines the Haar measure for the Unitary group $U(d)$.
Integration engine identifies variables via metadata tag :U.
IntegrateUnitary.dSU — Function
dSU(dim)Defines the Haar measure for the Special Unitary group $SU(d)$.
Integration for $SU(d)$ is performed via the $U(d)$ measure. For balanced polynomials, the integrals over $SU(d)$ and $U(d)$ coincide. In the current implementation, non-balanced expressions are handled by the same phase invariance rule as $U(d)$ and evaluate to zero.
Orthogonal & Symplectic
IntegrateUnitary.dO — Function
dO(dim)Defines the Haar measure for the real Orthogonal group $O(d)$ with dimension dim. Integration engine identifies variables via metadata tag :O.
IntegrateUnitary.dSp — Function
dSp(dim)Defines the Haar measure for the Symplectic group $Sp(d)$. The dimension dim must be even. Integration engine identifies variables via metadata tag :Sp.
Circular Ensembles
See Circular Ensembles for dCOE, dCSE, dCUE.
Gaussian Ensembles
IntegrateUnitary.dGUE — Function
dGUE(dim)Gaussian Unitary Ensemble (GUE) measure.
IntegrateUnitary.dGOE — Function
dGOE(dim)Gaussian Orthogonal Ensemble (GOE) measure.
IntegrateUnitary.dGSE — Function
dGSE(dim)Gaussian Symplectic Ensemble (GSE) measure. dim must be even.
IntegrateUnitary.dGinUE — Function
dGinUE(dim)Complex Ginibre Ensemble (GinUE) measure.
IntegrateUnitary.dGinOE — Function
dGinOE(dim)Real Ginibre Ensemble (GinOE) measure.
IntegrateUnitary.dGinSE — Function
dGinSE(dim)Quaternionic/Symplectic Ginibre Ensemble (GinSE) measure. dim must be even.
Pure States
IntegrateUnitary.dPsi — Function
dPsi(dim)Defines the Fubini-Study measure for a random pure state |psi> distributed according to the Haar measure.
Integration engine identifies variables via metadata tag :psi.
See Stiefel Manifolds for dStiefel.
Permutation Groups
IntegrateUnitary.dPerm — Function
dPerm(dim)Defines the Haar measure for the Symmetric group $S_d$ (permutation matrices) of dimension dim.
Integration engine identifies variables via metadata tag :Perm.
IntegrateUnitary.dCPerm — Function
dCPerm(dim)Defines the measure for Centered Permutation matrices $Y = P - J/d$ where $P \in S_d$.
Diagonal Unitary Matrices
IntegrateUnitary.dDiagUnitary — Function
dDiagUnitary(dim)Defines the measure for the group of diagonal unitary matrices (the torus $T^d$) of dimension dim. Integration engine identifies variables via metadata tag :DiagUnitary.
IntegrateUnitary.DiagonalUnitaryMeasure — Type
DiagonalUnitaryMeasure{D}Represents integration over the group of diagonal unitary matrices (the torus $T^d$). For a diagonal unitary matrix $V$, the only non-zero entries are $V_{ii} = e^{i\theta_i}$. Integration over $T^d$ is equivalent to independent phase integrations for each diagonal entry.
Symbolic Helpers
IntegrateUnitary.SymbolicMatrix — Type
SymbolicMatrix(name::Symbol)
SymbolicMatrix(name::Symbol, special_type::Symbol)
SymbolicMatrix(name::Symbol, special_type::Symbol, dim)
SymbolicMatrix(name::Symbol, is_adj::Bool, special_type::Symbol, dim)A wrapper associated with a symbolic name to represent a matrix in a coordinate-free way. Used in the symbolic trace logic (via tr_lazy) and for metadata-driven element-wise integration.
LinearAlgebra.tr — Function
tr(A::SymbolicMatrix)
tr(A::SymbolicMatrixProduct)Symbolic trace of a coordinate-free matrix expression. Returns a LazyTrace object that can be integrated.
IntegrateUnitary.symbolic_unitary — Function
symbolic_unitary(name::Symbol, d) -> SymbolicMatrixConstruct a d×d symbolic Haar-random unitary matrix named name. The returned object behaves as an AbstractMatrix and can be passed directly to integrate or used inside @integrate. d may be an integer or a Symbolics.Num symbolic variable for fully symbolic results.
Examples
using IntegrateUnitary, Symbolics
@variables d
U = symbolic_unitary(:U, d)
integrate(abs(U[1,1])^2, dU(d)) # 1/dIntegrateUnitary.symbolic_orthogonal — Function
symbolic_orthogonal(name::Symbol, d) -> SymbolicMatrixConstruct a d×d symbolic Haar-random orthogonal matrix named name. Use with dO(d) for integration over the orthogonal group O(d). d may be an integer or a Symbolics.Num symbolic variable.
Examples
using IntegrateUnitary, Symbolics
@variables d
O = symbolic_orthogonal(:O, d)
integrate(O[1,1]^2, dO(d)) # 1/dIntegrateUnitary.symbolic_symplectic — Function
symbolic_symplectic(name::Symbol, d) -> SymbolicMatrixConstruct a d×d symbolic Haar-random symplectic matrix named name. Use with dSp(d) for integration over the symplectic group Sp(d). d must be even (validated at integration time).
Examples
using IntegrateUnitary, Symbolics
Sp = symbolic_symplectic(:Sp, 4)
integrate(abs(Sp[1,1])^2, dSp(4)) # 1/4IntegrateUnitary.symbolic_pure_state — Function
symbolic_pure_state(name::Symbol, d) -> SymbolicMatrixConstruct a symbolic Haar-random pure state (column vector) in ℂ^d, named name. The returned object has dimensions (d, 1) and is indexed as psi[i, 1]. Use with dPsi(d) for integration over the Fubini–Study measure.
Examples
using IntegrateUnitary, Symbolics
@variables d
psi = symbolic_pure_state(:psi, d)
integrate(abs(psi[1,1])^2, dPsi(d)) # 1/dIntegrateUnitary.symbolic_permutation — Function
symbolic_permutation(name::Symbol, d) -> SymbolicMatrixConstruct a symbolic random d×d permutation matrix named name. Use with dPerm(d) for integration over the uniform measure on the symmetric group S_d.
Examples
using IntegrateUnitary, Symbolics
@variables d
P = symbolic_permutation(:P, d)
integrate(P[1,1] * P[2,2], dPerm(d)) # 1/(d*(d-1))Quantum Information Utilities
IntegrateUnitary.partial_trace — Function
partial_trace(M, dims, subsystem)Compute the partial trace of matrix M over the specified subsystem. dims is a tuple/vector of dimensions for each subsystem. subsystem is the index of the subsystem to be TRACED OUT.
Example: For a bipartite system with dims=(2, 2), partial_trace(M, (2, 2), 2) returns the reduced density matrix of the first subsystem.
Internal / Advanced
These functions are part of the internal machinery but documented for development reference.
Integration Engine & Helpers
IntegrateUnitary.integrate_indices — Function
integrate_indices(U_idxs, U_bar_idxs, dim)Low-level integration function using Weingarten calculus (Unitary).
IntegrateUnitary.integrate_indices_symplectic — Function
integrate_indices_symplectic(indices, dim)Low-level integration function using Symplectic Weingarten calculus. Formula: sum{pi, sigma} Jpi(i) * J_sigma(j) * Wg^Sp(pi, sigma)
IntegrateUnitary.integrate_indices_permutation — Function
integrate_indices_permutation(indices, dim)Integration over the symmetric group S_d.
IntegrateUnitary.integrate_indices_diagonal — Function
integrate_indices_diagonal(U_idxs, U_bar_idxs, dim)Low-level integration function for the Diagonal Unitary group (Torus). The integral is non-zero (equal to 1) only if the multisets of indices of U and Ubar are identical. Indices are already checked to be diagonal (i == j) in processterm.
IntegrateUnitary.integrate_indices_coe — Function
integrate_indices_coe(indices, U_bar_indices, dim)Integration of COE terms by reducing to Haar integration.
Index Mapping
$S = U U^T$, so each $S_{ij}$ expands as:
\[S_{ij} = \sum_k U_{ik} U_{jk}, \qquad \bar{S}_{ij} = \sum_k \bar{U}_{ik} \bar{U}_{jk}\]
Each $S_{i_k j_k}$ produces two U-type indices $(i_k, a_k)$ and $(j_k, a_k)$ with shared dummy column $a_k$. Similarly each $\bar{S}_{p_k q_k}$ produces two $\bar{U}$-type indices with shared dummy column $b_k$.
Weingarten Integration
With $m$ S-terms and $m$ $\bar{S}$-terms, we have $2m$ U and $2m$ $\bar{U}$ indices. The Haar integral gives:
\[\sum_{\sigma, \tau \in S_{2m}} \delta_{\text{rows}}(\sigma) \cdot d^{\#\text{loops}(\tau)} \cdot \text{Wg}(\sigma\tau^{-1}, d)\]
The dummy column summation $\sum_{a,b} \delta_{\text{cols}}(\tau)$ yields $d^{\#\text{loops}(\tau)}$, where the number of loops is the number of connected components in a bipartite graph:
- Vertices: $a_1 \ldots a_m$ and $b_1 \ldots b_m$
- Structural edges: $a_k$ pairs columns $(2k{-}1, 2k)$; same for $b_k$
- $\tau$-edges: $a_{\lceil r/2 \rceil} \leftrightarrow b_{\lceil \tau(r)/2 \rceil}$
The loop count is computed via union-find on these $2m$ variable nodes.
IntegrateUnitary.integrate_indices_gue — Function
integrate_indices_gue(indices, dim)Low-level integration function using Wick's theorem for GUE. Formula: sum{pi in PairPartitions} prod{(u, v) in pi} delta(iu, jv) * delta(ju, iv)
IntegrateUnitary.integrate_indices_goe — Function
integrate_indices_goe(indices, dim)Low-level integration function using Wick's theorem for GOE. Formula: sum{pi in PairPartitions} prod{(u, v) in pi} (delta(iu, kv)delta(ju, lv) + delta(iu, lv)delta(ju, kv)) where pair is H{iu ju} and H{kv lv} (indices re-labeled for clarity).
IntegrateUnitary.integrate_indices_ginue — Function
integrate_indices_ginue(u_indices, ub_indices, dim)Low-level integration for Complex Ginibre Ensemble. Formula: sum{sigma in Sn} prodm delta(im, ibsigma(m)) * delta(jm, jb_sigma(m))
IntegrateUnitary.integrate_indices_ginoe — Function
integrate_indices_ginoe(indices, dim)Low-level integration for Real Ginibre Ensemble. Formula: sum{pi in PairPartitions} prod{(u, v) in pi} delta(iu, iv) * delta(ju, jv)
IntegrateUnitary.integrate_indices_ginse — Function
integrate_indices_ginse(indices, dim)Low-level integration for Symplectic Ginibre Ensemble. Uses duality with GinOE.
IntegrateUnitary.tr_lazy — Function
tr_lazy(product) -> LazyTraceBuild a LazyTrace representing the symbolic trace of a matrix product without evaluating any indices. product may be a SymbolicMatrix, a SymbolicMatrixProduct, or an AbstractVector of matrix factors.
Prefer the tr() overload when constructing expressions with * (it delegates here automatically). Call tr_lazy directly when you already hold a pre-assembled Vector of matrix factors.
Examples
using IntegrateUnitary, Symbolics
@variables d
U = symbolic_unitary(:U, d)
A = SymbolicMatrix(:A)
# Equivalent ways to build tr(U*A*U'):
expr1 = tr(U * A * U')
expr2 = tr_lazy([U, A, U'])
integrate(expr1, dU(d)) # tr(A)/dSee also: tr, LazyTrace, symbolic_unitary
IntegrateUnitary.LazyTrace — Type
LazyTrace(cycles::Vector{Vector{AbstractMatrix}}, prefactor::Union{Num, Number})A lazy representation of a trace (or product of traces) of matrix products. Used to represent expressions like tr(A*B) * tr(C) symbolically before integration.
IntegrateUnitary.LazySum — Type
LazySum(terms::Vector{LazyTrace})A lazy representation of a sum of LazyTrace objects. Enables symbolic integration of expressions like tr(A*B) + tr(C*D).
IntegrateUnitary.check_library — Function
check_library(expr, measure)Check if the integral of expr over measure is available in the pre-computed library. Returns the symbolic result if found, otherwise nothing.
IntegrateUnitary.tr_val — Function
tr_val(factors::AbstractVector)Evaluates the trace of a product of matrices. If all factors are concrete, returns the numeric trace. If any factor is symbolic, returns a symbolic representation normalized by circular shifts and adjoints to ensure unique naming.
IntegrateUnitary._expand_asymptotic — Function
_expand_asymptotic(ex, d, order)Helper to expand a rational function of d in powers of 1/d. Uses a series dictionary representation for robust recursive processing.
IntegrateUnitary._poly_degree — Function
_poly_degree(p, d)Helper to get the degree of a polynomial p in variable d.
IntegrateUnitary._ensure_symbolic_dim — Function
_ensure_symbolic_dim(d)Ensure dimension d is a proper Symbolics variable. If d unwraps to a plain Symbol, wrap it via Symbolics.variable; otherwise return as-is.
IntegrateUnitary._try_numeric — Function
_try_numeric(v)Attempt to convert a value to a clean numeric form. Returns the converted value, or nothing if conversion is not possible.
AbstractFloat→ rationalizedReal→ returned as-is
IntegrateUnitary._try_extract_int — Function
_try_extract_int(x)Extract a plain Int from x if it wraps a concrete integer value. Returns the Int, or nothing if x is symbolic or non-integer. Handles: plain Integer, Num-wrapped integer constants (e.g. Num(2)).
IntegrateUnitary.robust_substitute — Function
robust_substitute(ex, dict)Substitute symbolic variables in ex using dict. Falls back to manual Postwalk traversal with unwrapped keys if Symbolics.substitute doesn't produce a different expression.
IntegrateUnitary.get_full_cycle_type — Function
get_full_cycle_type(pi, sigma)Returns the cycle type of the union of two pair partitions as a sorted partition of k. The union forms cycles of lengths 2l1, 2l2, ... where sum li = k. Returns [l1, l_2, ...] sorted descending.
IntegrateUnitary.get_weingarten_reduced_data — Function
get_weingarten_orthogonal_data(k, d)Internal function to generate the Orthogonal Weingarten matrix. The matrix $G$ is a Gram matrix of size $(2k-1)!! \times (2k-1)!!$ where $G_{\pi, \sigma} = d^{\ell(\pi, \sigma)}$. The Weingarten matrix is the inverse of $G$.
Matcher and Logic
IntegrateUnitary.AbstractIndexMatcher — Type
AbstractIndexMatcherAbstract base type for strategies that identify random-matrix entries inside a symbolic expression. Concrete subtypes must implement
match_index(m::AbstractIndexMatcher, t) -> Union{Tuple, Nothing}returning (tag, i, j) when t is a recognised random-matrix entry (where tag is a Symbol like :U or :U_bar, and i, j are row/column indices or nothing for metadata-only matching), or nothing otherwise.
The built-in subtype is MetadataMatcher.
IntegrateUnitary.MetadataMatcher — Type
MetadataMatcher(type_tag::Symbol)An AbstractIndexMatcher that recognises random-matrix entries by the special_type metadata attached to a SymbolicMatrix symbol.
type_tag must match the special_type field of the target matrix, e.g.:
:Ufor Haar-unitary matrices (symbolic_unitary):Ofor orthogonal matrices (symbolic_orthogonal):Spfor symplectic matrices (symbolic_symplectic):psifor pure-state vectors (symbolic_pure_state)
Conjugate entries (is_adj = true) are tagged as Symbol(type_tag, :_bar) (e.g. :U_bar).
IntegrateUnitary._integrate_core — Function
_integrate_core(expr, dim, subs_dict, matcher, measure_type=:U)The internal integration engine. It performs several steps:
- Normalization/Rewriting: Expands
abs2(z),real(z),imag(z)into explicit polynomials. - Substitution: Replaces symbolic variables with internal atomic representatives via
subs_dict. - Expansion: Distributes products over sums to get a sum of monomials.
- Monomial Integration: For each monomial, invokes
process_termto identify indices and apply the appropriate integration rule (Weingarten or Wick).
IntegrateUnitary.process_term — Function
process_term(term, matcher, dim, measure_type)Integrates a single monomial term.
- Index Collection: Traverses the term to find all random matrix elements using the
matcher. - Dispatch: Calls specific index integration functions based on
measure_type::U:integrate_indices(Weingarten):O:integrate_indices_orthogonal:Sp:integrate_indices_symplectic:GUE,:GOE,:GSE: Wick contraction rules.
IntegrateUnitary.weingarten — Function
weingarten(partition_type::Vector{Int}, d)Computes the Unitary Weingarten function \text{Wg}(\sigma, d) where \sigma is a permutation with cycle type given by partition_type.
The Weingarten function is defined as the sum over irreducible representations of S_n:
\[\text{Wg}(\sigma, d) = \frac{1}{(n!)^2} \sum_{\lambda \vdash n, \ell(\lambda) \le d} \frac{(f^\lambda)^2 \chi^\lambda(\sigma)}{s_\lambda(1^d)}\]
where f^\lambda is the dimension of the Sn irrep, \chi^\lambda(\sigma) is the character, and s\lambda(1^d) is the dimension of the U(d) irrep.
Reference:
- Collins, B., & Śniady, P. (2006). Integration with respect to the Haar measure on unitary, orthogonal and symplectic groups. Communications in Mathematical Physics.
IntegrateUnitary.ParityUnionFind — Type
ParityUnionFindUnion-find data structure with parity tracking for CSE integration. Tracks whether the path from a node to its root has odd or even parity.
Weingarten & Combinatorics
IntegrateUnitary.weingarten_orthogonal_val — Function
weingarten_orthogonal_val(pi, sigma, d)Returns the Orthogonal Weingarten function value \text{Wg}^O(\pi, \sigma, d).
Reference:
- Collins, B., & Śniady, P. (2006). Integration with respect to the Haar measure on unitary, orthogonal and symplectic groups.
IntegrateUnitary.weingarten_symplectic_val — Function
weingarten_symplectic_val(pi, sigma, d)Returns the Symplectic Weingarten function value \text{Wg}^{Sp}(\pi, \sigma, d). Uses the duality relation:
\[\text{Wg}^{Sp}(\pi, \sigma, d) = (-1)^{\text{loops}(\pi, \sigma)} \text{Wg}^{O}(\pi, \sigma, -d)\]
where loops is the number of cycles in the union of the two pair partitions.
Reference:
- Collins, B., & Śniady, P. (2006). Integration with respect to the Haar measure on unitary, orthogonal and symplectic groups.
IntegrateUnitary.get_pair_partitions — Function
get_pair_partitions(n)Generate all partitions of the set {1, ..., n} into n/2 disjoint pairs. n must be even. The number of such partitions is given by the double factorial (n-1)!!.
These partitions are also known as perfect matchings of the complete graph K_n.
IntegrateUnitary.get_matching_pair_partitions_filtered — Function
get_matching_pair_partitions_filtered(indices)Returns all pair partitions of 1..2k such that for every pair (a,b), indices[a] == indices[b].
IntegrateUnitary.canonicalize_pair_partition — Function
canonicalize_pair_partition(p)Returns a sorted list of pairs (min, max), sorted by min, to uniquely identify a matching.
IntegrateUnitary.conjugate_partition — Function
conjugate_partition(part::Vector{Int})Returns the conjugate partition \lambda' of a partition \lambda. The conjugate partition is obtained by transposing the Young diagram of \lambda. Mathematically, \lambda'i = \text{card}{j : \lambdaj \ge i}.
IntegrateUnitary.murnaghan_nakayama — Function
murnaghan_nakayama(lambda, mu)Computes the character \chi^\lambda(\mu) of the symmetric group S_n for the irrep \lambda at the conjugacy class with cycle type \mu using the Murnaghan-Nakayama rule.
The rule states:
\[\chi^\lambda(\mu) = \sum_{T \in RIM(\lambda, \mu)} (-1)^{\text{ht}(T)}\]
where the sum is over "rim hook" tableaux of shape \lambda and content \mu.
Reference:
- Murnaghan, F. D. (1937). The characters of the symmetric group.
- Nakayama, T. (1940). On some finite group of substitutions.
IntegrateUnitary.character_at_id — Function
character_at_id(part::Vector{Int})Calculates the character of the symmetric group S_n at the identity element for the irreducible representation (irrep) corresponding to the partition part (\lambda). This is equivalent to the dimension f^\lambda of the irrep.
The dimension is given by the Hook Length Formula:
\[f^\\lambda = \\frac{n!}{\\prod_{(i,j) \\in \\lambda} h_{\\lambda}(i,j)}\]
where h_\lambda(i,j) is the hook length of the cell (i,j) in the Young diagram of \lambda.
Reference:
- Frame, J. S., Robinson, G. de B., & Thrall, R. M. (1954). The hook graphs of the symmetric group.
IntegrateUnitary.irrep_dimension — Function
irrep_dimension(part::Vector{Int}, d)Computes the dimension s_\lambda(1^d) of the irreducible representation of the unitary group U(d) (or the Schur polynomial at ones) corresponding to the partition \lambda.
The dimension is given by the Hook-Content Formula:
\[\text{dim}_d(\lambda) = \prod_{(i,j) \in \lambda} \frac{d + j - i}{h_\lambda(i,j)}\]
where h_\lambda(i,j) is the hook length and j-i is the content of the cell (i,j).
This implementation supports symbolic dimension d, returning a rational function in d.
Reference:
- Stanley, R. P. (1999). Enumerative Combinatorics, Vol. 2.
IntegrateUnitary.compute_symplectic_contraction — Function
compute_symplectic_contraction(partition, indices, dim)Computes prod_{(u, v) in partition} J(indices[u], indices[v]). J = [0 I; -I 0].
IntegrateUnitary.weingarten_orthogonal_val_canonical — Function
weingarten_orthogonal_val_canonical(c_pi, c_sigma, d)Internal version of weingarten_orthogonal_val that assumes arguments are already canonical.
IntegrateUnitary.INTEGRATION_RULES — Constant
INTEGRATION_RULESA dictionary mapping measure types (symbols) to their respective integration rule functions. Each rule function should have the signature (u_indices, u_bar_indices, dim, measure_type).