API Reference
This page provides a complete reference for all exported functions in FibonacciChain.jl.
Basis and Hamiltonian Functions
FibonacciChain.anyon_basis — Function
anyon_basis(model::AnyonModel; symmetry_block=nothing)
anyon_basis(anyon_type::AbstractAnyonType, ::Type{T}; pbc::Bool=true, symmetry_block=nothing) where {N, T <: BitStr{N}}Generate basis states for an anyon chain model.
This function provides two interfaces: a high-level one that takes an AnyonModel object, and a low-level one that takes the model parameters directly.
Arguments
High-level interface
model::AnyonModel: AnAnyonModelobject containing the system parameters (anyon type, size, and boundary conditions).symmetry_block: Optional. The topological charge sector to filter the basis. Not yet implemented.
Low-level interface
anyon_type::AbstractAnyonType: The type of anyon, e.g.,FibonacciAnyon()orIsingAnyon().T::Type: TheBitStrtype specifying the chain lengthN.pbc::Bool=true: Specifies whether to use periodic (true) or open (false) boundary conditions.symmetry_block: Optional. The topological charge sector.
Returns
Vector{T}: A sorted vector of the basis states.
Examples
julia> using FibonacciChain, BitBasis
julia> model = AnyonModel(FibonacciAnyon(), 4; pbc=true); basis = anyon_basis(model); length(basis)
7
julia> basis_fibo = anyon_basis(FibonacciAnyon(), BitStr{4, Int}; pbc=true); length(basis_fibo)
7
julia> model_ising = AnyonModel(IsingAnyon(), 4; pbc=true); basis_ising = anyon_basis(model_ising); length(basis_ising)
16anyon_basis(AT::AbstractAnyonType, ::Type{T}, k::Int; symmetry_block=nothing) where {N, T <: BitStr{N}}
anyon_basis(model::AnyonModel, k::Int; symmetry_block=nothing)Generate basis states in specific momentum sector k and topological sector symmetry_block.
Arguments
Low-level interface
AT::AbstractAnyonType: Anyon type instance (e.g.,FibonacciAnyon(),IsingAnyon())T::Type: BitStr type specifying chain length Nk::Int: Momentum sector (0 ≤ k ≤ N-1)symmetry_block: Topological charge sector (optional)
High-level interface
model::AnyonModel: Anyon model containing system parametersk::Int: Momentum sector (0 ≤ k ≤ N-1)symmetry_block: Topological charge sector (optional)
Returns
Vector{T}: Basis states in momentum sector kDict{T, Vector{T}}: Representative mapping for translation equivalence classes
Examples
julia> using FibonacciChain, BitBasis
julia> model = AnyonModel(FibonacciAnyon(), 6; pbc=true);
julia> basisK, basis_dic = anyon_basis(model, 0);
julia> length(basisK) > 0
trueFibonacciChain.anyon_ham — Function
anyon_ham(model::AnyonModel; kwargs...)
anyon_ham(model::AnyonModel, ::Type{T}; kwargs...) where {N, T <: BitStr{N}}Construct the Hamiltonian matrix for a 1D anyon chain.
Arguments
High-level interface
model::AnyonModel: AnAnyonModelobject specifying the anyon type, system size, boundary conditions, and interaction type.kwargs...: Additional model-specific parameters, such asJandhfor the Ising model.
Low-level usage (discouraged)
While not a direct low-level function, you can construct the model on the fly: anyon_ham(AnyonModel(IsingAnyon(), N), ...)
Returns
Matrix{Float64}: The Hamiltonian matrix constructed in the chosen basis.
Supported Models
- Fibonacci Anyons: Supports
:Antiferroand:Ferrointeraction terms. - Ising Anyons: Transverse field Ising model.
Examples
julia> using FibonacciChain, BitBasis
julia> N = 4; model_fibo = AnyonModel(FibonacciAnyon(), N; pbc=true, measure_operator=:Antiferro);
julia> H_fibo = anyon_ham(model_fibo); size(H_fibo)
(7, 7)
julia> model_ising = AnyonModel(IsingAnyon(), N; pbc=true); H_ising = anyon_ham(model_ising; J=1.0, h=1.0); size(H_ising)
(16, 16)anyon_ham(model::AnyonModel, ::Type{T}, k::Int; symmetry_block=nothing, kwargs...) where {N, T <: BitStr{N}}
anyon_ham(model::AnyonModel, k::Int; symmetry_block=nothing, kwargs...)Construct Hamiltonian matrix in specific momentum sector for 1D anyon chain.
Arguments
model::AnyonModel: Anyon model containing system parametersT::Type: BitStr type specifying chain length N (optional)k::Int: Momentum sector (0 ≤ k ≤ N-1)symmetry_block: Topological charge sector (optional)kwargs...: Additional model parameters, e.g.,J,hfor Ising model
Returns
Matrix: Hamiltonian matrix in chosen momentum sector (ComplexF64, or real if k=0 or k=N/2)
Examples
julia> using FibonacciChain
julia> model = AnyonModel(FibonacciAnyon(), 6; pbc=true);
julia> H_k0 = anyon_ham(model, 0);
julia> size(H_k0)[1] > 0
trueanyon_ham(model::AnyonModel, sites::Vector{<:Index}; kwargs...)Construct anyon chain Hamiltonian as Matrix Product Operator (MPO).
Arguments
model::AnyonModel: Anyon model containing system parameterssites::Vector{<:Index}: ITensor site indiceskwargs...: Additional parameters (e.g.,J,hfor Ising model)
Returns
MPO: Hamiltonian as Matrix Product Operator- Fibonacci: three-body interactions based on Fibonacci fusion rules
- Ising: transverse field Ising model with ZZ and X terms
FibonacciChain.anyon_ham_sparse — Function
anyon_ham_sparse(model::AnyonModel; kwargs...)Construct the anyon chain Hamiltonian as a sparse matrix.
Arguments
model::AnyonModel: AnAnyonModelobject containing the system parameters.kwargs...: Additional keyword arguments passed toactingHam, e.g.,J,hfor the Ising model.
Returns
SparseMatrixCSC{Float64, Int}: The sparse Hamiltonian matrix in the anyon basis.
Examples
julia> using FibonacciChain, BitBasis, SparseArrays, LinearAlgebra
julia> N = 6; model = AnyonModel(FibonacciAnyon(), N; pbc=true); H_sparse = anyon_ham_sparse(model);
julia> H_sparse isa SparseMatrixCSC
true
julia> ishermitian(H_sparse)
true
julia> H_dense = anyon_ham(model); norm(Matrix(H_sparse) - H_dense) < 1e-10
trueanyon_ham_sparse(model::AnyonModel, k::Int; symmetric_block=nothing, kwargs...)Construct the Hamiltonian in a specific momentum and/or topological sector as a sparse matrix.
Arguments
model::AnyonModel: AnAnyonModelobject containing the system parameters.k::Int: The momentum quantum number (0 ≤ k ≤ N-1).symmetric_block: Optional. The topological charge sector (e.g.,:trivial,:nontrivial).kwargs...: Additional keyword arguments passed toactingHam, e.g.,J,h.
Returns
SparseMatrixCSC: The sparse Hamiltonian in the specified symmetric sector. The element type will beComplexF64for generalkandFloat64fork=0ork=N/2.
FibonacciChain.mapst_sec2tot — Function
mapst_sec2tot(AT::AbstractAnyonType, ::Type{T}, state::Vector{ET}, k::Int; pbc::Bool=true) where {N, T <: BitStr{N}, ET}
mapst_sec2tot(model::AnyonModel, state::Vector{ET}, k::Int)Map state in momentum sector to total Hilbert space.
Arguments
High-level interface
model::AnyonModel: AnAnyonModelobject containing system parametersstate::Vector{ET}: The state vector in the Hilbert space of the specified momentum sectork::Int: The momentum sector index (0 ≤ k ≤ N-1)
Low-level interface
AT::AbstractAnyonType: The anyon type.T::Type: TheBitStrtype specifying the chain lengthN.state::Vector{ET}: The state vector in the momentum sector.k::Int: The momentum sector index.pbc::Bool=true: Specifies periodic boundary conditions (required for momentum sectors).
Returns
Vector{ET}: The state vector represented in the total Hilbert space basis.
Examples
julia> using FibonacciChain, LinearAlgebra
julia> model = AnyonModel(FibonacciAnyon(), 6; pbc=true); H_k0 = anyon_ham(model, 0);
julia> eigvecs_k0 = eigvecs(H_k0); kstate = eigvecs_k0[:, 1];
julia> total_state = mapst_sec2tot(model, kstate, 0); length(total_state) == length(anyon_basis(model))
trueReduced Density Matrices
FibonacciChain.anyon_rdm — Function
anyon_rdm(AT::AbstractAnyonType, ::Type{T}, subsystems::Vector{Int64}, state::Vector{ET}; pbc::Bool=true) where {N, T <: BitStr{N}, ET}
anyon_rdm(AT::AbstractAnyonType, ::Type{T}, subsystems::Vector{Int64}, state::Matrix{ET}; pbc::Bool=true) where {N, T <: BitStr{N}, ET}
anyon_rdm(model::AnyonModel, subsystems::Vector{Int64}, state::Union{Vector{ET}, Matrix{ET}})Compute the reduced density matrix (RDM) for a specified subsystem.
This function can operate on both state vectors and density matrices.
Arguments
High-level interface
model::AnyonModel: AnAnyonModelobject containing the system parameters.subsystems::Vector{Int64}: A vector of site indices to keep in the reduced system. Indices are 1-based.state::Union{Vector{ET}, Matrix{ET}}: The quantum state, either as a state vector or a full density matrix.
Low-level interface
AT::AbstractAnyonType: The anyon type, e.g.,FibonacciAnyon().T::Type: TheBitStrtype specifying the total chain lengthN.subsystems::Vector{Int64}: The indices of the subsystem sites.state::Union{Vector{ET}, Matrix{ET}}: The quantum state vector or density matrix.pbc::Bool=true: Specifies whether the original system has periodic boundary conditions.
Returns
Matrix{ET}: The reduced density matrix for the specified subsystem.
Details
The subsystem indices are 1-based and refer to the site positions in the chain. The function correctly handles both contiguous and non-contiguous subsystems.
Examples
julia> using FibonacciChain, BitBasis, LinearAlgebra, Random
julia> Random.seed!(42); N = 4; model = AnyonModel(FibonacciAnyon(), N; pbc=true);
julia> basis = anyon_basis(model); state = randn(ComplexF64, length(basis)); state ./= norm(state);
julia> rdm = anyon_rdm(model, [1, 2], state); size(rdm)
(3, 3)
julia> ishermitian(rdm)
true
julia> isapprox(tr(rdm), 1.0; atol=1e-10)
trueFibonacciChain.anyon_rdm_sec — Function
anyon_rdm_sec(model::AnyonModel, subsystems::Vector{Int64}, kstate::Vector{ET}, k::Int) where {ET}Compute reduced density matrix for specified subsystems from quantum state in specific momentum sector.
Arguments
model::AnyonModel: Anyon model containing system parameterssubsystems::Vector{Int64}: Indices of subsystem sites to keepkstate::Vector{ET}: State vector in momentum sector Hilbert spacek::Int: Momentum sector (0 ≤ k ≤ N-1)
Returns
Matrix{ET}: Reduced density matrix in total Hilbert basis
Examples
julia> using FibonacciChain, LinearAlgebra
julia> model = AnyonModel(FibonacciAnyon(), 6; pbc=true);
julia> H_k0 = anyon_ham(model, 0);
julia> eigvecs_k0 = eigvecs(H_k0);
julia> kstate = eigvecs_k0[:, 1];
julia> rdm = anyon_rdm_sec(model, [1, 2, 3], kstate, 0);
julia> isapprox(tr(rdm), 1.0; atol=1e-10)
trueFibonacciChain.disjoint_rdm — Function
disjoint_rdm(AT1::AbstractAnyonType, AT2::AbstractAnyonType, ::Type{T1}, ::Type{T2}, subsystemsA::Vector{Int64}, subsystemsB::Vector{Int64}, state::Vector{ET}; totalsubApbc::Bool=false, totalsubBpbc::Bool=false) where {N1, N2, T1 <: BitStr{N1}, T2 <: BitStr{N2}, ET}
disjoint_rdm(model1::AnyonModel, model2::AnyonModel, subsystemsA::Vector{Int64}, subsystemsB::Vector{Int64}, state::Vector{ET})Compute reduced density matrix for two joint different systems, or two parallel chains.
Arguments
Low-level interface
AT1::AbstractAnyonType: Anyon type for subsystem AAT2::AbstractAnyonType: Anyon type for subsystem BT1::Type: BitStr type specifying chain length N1T2::Type: BitStr type specifying chain length N2subsystemsA::Vector{Int64}: Indices of subsystem A sites to keepsubsystemsB::Vector{Int64}: Indices of subsystem B sites to keepstate::Vector{ET}: Quantum state vector in total Hilbert spacetotalsubApbc::Bool=false: Whether subsystem A is periodictotalsubBpbc::Bool=false: Whether subsystem B is periodic
High-level interface
model1::AnyonModel: The model for the first subsystem (A).model2::AnyonModel: The model for the second subsystem (B).subsystemsA::Vector{Int64}: Indices of subsystem A sites to keepsubsystemsB::Vector{Int64}: Indices of subsystem B sites to keepstate::Vector{ET}: Quantum state vector in total Hilbert space
Returns
Matrix{ET}: Reduced density matrix
Examples
julia> using FibonacciChain, BitBasis, LinearAlgebra
julia> N1, N2 = 4, 4;
julia> model1 = AnyonModel(FibonacciAnyon(), N1; pbc=true);
julia> model2 = AnyonModel(FibonacciAnyon(), N2; pbc=true);
julia> basisA = anyon_basis(model1);
julia> basisB = anyon_basis(model2);
julia> state = randn(ComplexF64, length(basisA) * length(basisB));
julia> state ./= norm(state);
julia> rdm = disjoint_rdm(model1, model2, [1,2], [1,2], state);
julia> size(rdm) # reduced density matrix dimension
(9, 9)
julia> ishermitian(rdm) # should be Hermitian
true
julia> isapprox(tr(rdm), 1.0; atol=1e-10) # trace should be 1
trueEntanglement and Observables
FibonacciChain.ee — Function
ee(subrm::Matrix{ET}) where {ET}Calculate entanglement entropy from reduced density matrix.
Arguments
subrm::Matrix{ET}: Hermitian reduced density matrix
Returns
Float64: von Neumann entanglement entropy
Examples
julia> using FibonacciChain, LinearAlgebra
julia> # Create a simple 2x2 density matrix
ρ = [0.5 0.0; 0.0 0.5]; # Maximally mixed state
julia> entropy = ee(ρ);
julia> abs(entropy - log(2)) < 1e-10 # Should equal log(2) ≈ 0.693
true
julia> # Pure state has zero entropy
ρ_pure = [1.0 0.0; 0.0 0.0];
julia> ee(ρ_pure) ≈ 0.0
trueFibonacciChain.ee_mps — Function
ee_mps(ψ::MPS, b::Int)Calculate entanglement entropy of MPS state with bipartition at bond b.
Arguments
ψ::MPS: MPS stateb::Int: Bond index for bipartition (1 ≤ b < N)
Returns
Float64: Entanglement entropy (von Neumann entropy) at bond b
FibonacciChain.anyon_eelis — Function
anyon_eelis(model::AnyonModel, state::Union{Vector{ET}, Matrix{ET}}) where {ET}Calculate entanglement entropy profile along the chain for given quantum state or density matrix.
Arguments
model::AnyonModel: Anyon model containing system parameters (N, pbc, anyon_type)state::Union{Vector{ET}, Matrix{ET}}: Quantum state vector or density matrix in anyon basis
Returns
Vector{Float64}: Entanglement entropy at each bipartition from left to right
Examples
julia> using FibonacciChain, LinearAlgebra
julia> N = 6;
julia> model = AnyonModel(FibonacciAnyon(), N; pbc=true);
julia> # Create ground state of Fibonacci Hamiltonian
H = anyon_ham(model);
julia> eigenvals, eigenvecs = eigen(H);
julia> ground_state = eigenvecs[:, 1];
julia> # Calculate entanglement entropy profile
ee_profile = anyon_eelis(model, ground_state);
julia> length(ee_profile) == N - 1 # Profile has N-1 points
true
julia> all(x -> x ≥ 0, ee_profile) # All entropies are non-negative
trueanyon_eelis(model::AnyonModel, ψ::MPS)Calculate entanglement entropy profile along the chain for MPS state.
Arguments
model::AnyonModel: Anyon model containing system parametersψ::MPS: MPS state
Returns
Vector{Float64}: Entanglement entropy at each bipartition from left to right
Examples
julia> using FibonacciChain, ITensorMPS, ITensors
julia> model = AnyonModel(FibonacciAnyon(), 6; pbc=true);
julia> ψ_gs, E0 = anyon_mps_gst(model, maxdim=10, outputlevel=0);
julia> ee_profile = anyon_eelis(model, ψ_gs);
julia> length(ee_profile) == model.N - 1 # Profile has N-1 points
true
julia> all(x -> x ≥ 0, ee_profile) # All entropies are non-negative
trueSymmetry and Topological Operations
FibonacciChain.translation_matrix — Function
translation_matrix(model::AnyonModel)Generate translation operator matrix for anyon basis states.
Arguments
model::AnyonModel: Anyon model containing system parameters (N, pbc, anyon_type)
Returns
Matrix{Float64}: Translation matrix mapping each basis state to its translated version
FibonacciChain.inversion_matrix — Function
inversion_matrix(model::AnyonModel)Generate spatial inversion operator matrix for anyon basis states.
Arguments
model::AnyonModel: Anyon model containing system parameters (N, pbc, anyon_type)
Returns
Matrix{Float64}: Inversion matrix mapping each basis state to its spatially reflected version
FibonacciChain.braidingsqmap — Function
braidingsqmap(model::AnyonModel{FibonacciAnyon}, state::Vector{ET}, idx::Int) where {ET}Apply braiding squared operation to quantum state at specified site.
Arguments
model::AnyonModel{FibonacciAnyon}: Fibonacci anyon model containing system parameters (N, pbc)state::Vector{ET}: Quantum state vector in anyon basisidx::Int: Site index for braiding operation
Returns
Vector{ComplexF64}: Transformed state after braiding operation
FibonacciChain.topological_symmetry_basismap — Function
topological_symmetry_basismap(model::AnyonModel{FibonacciAnyon}, state::T) where {N, T <: BitStr{N}}Compute topological symmetry coefficients for all basis states relative to given state.
Arguments
model::AnyonModel{FibonacciAnyon}: Fibonacci anyon modelstate::T: Reference state configuration
Returns
Vector{Float64}: Coefficients for each basis state
Examples
julia> using FibonacciChain, BitBasis
julia> N = 4; T = BitStr{N, Int};
julia> model = AnyonModel(FibonacciAnyon(), N; pbc=true);
julia> state = T(0b0000); # Vacuum state
julia> coeffs = topological_symmetry_basismap(model, state);
julia> length(coeffs) == length(anyon_basis(model))
trueFibonacciChain.Fsymmetry_coef — Function
Fsymmetry_coef(::FibonacciAnyon, ::Type{T}, state::T, base::T; pbc::Bool=true) where {N, T <: BitStr{N}}
Fsymmetry_coef(model::AnyonModel{FibonacciAnyon}, state::T, base::T)Compute topological symmetry coefficient for state in given base configurations for Fibonacci anyon chain.
Arguments
Low-level interface
::FibonacciAnyon: Fibonacci anyon type instance::Type{T}: BitStr type specifying chain length Nstate::T: Target state configurationbase::T: Base state configurationpbc::Bool=true: Periodic boundary conditions
High-level interface
model::AnyonModel{FibonacciAnyon}: Anyon model containing system parametersstate::T: Target state configurationbase::T: Base state configuration
Returns
Float64: Topological symmetry coefficient based on Fibonacci fusion rules
Examples
julia> using FibonacciChain, BitBasis
julia> N = 4; T = BitStr{N, Int};
julia> model = AnyonModel(FibonacciAnyon(), N; pbc=true);
julia> state = T(0b1010); ϕ = (1 + sqrt(5)) / 2;
julia> base = T(0b0101);
julia> coef = Fsymmetry_coef(model, state, base); coef ≈ 0.3819660112501051
trueLadder System Functions
FibonacciChain.ladderChoi — Function
ladderChoi(model::AnyonModel{FibonacciAnyon}, p::Float64, state::Vector{ET}) where {ET}Apply probabilistic braiding noise channel to density matrix state in vec form.
Arguments
model::AnyonModel{FibonacciAnyon}: Fibonacci anyon model containing system parameters (N, pbc)p::Float64: Braiding probability (0 ≤ p ≤ 1)state::Vector{ET}: Density matrix state in vectorized form
Returns
Vector{ET}: State after applying Choi map with braiding noise
Implements noise channel: ρ → (1-p)ρ + p B(ρ) where B is braiding operation.
FibonacciChain.ladderrdm — Function
ladderrdm(model::AnyonModel, subsystems::Vector{Int64}, state::Vector{ET}) where {ET}Compute reduced density matrix for specified subsystem from vectorized density matrix.
Arguments
model::AnyonModel: Anyon model containing system parameters (N, pbc, anyon_type)subsystems::Vector{Int64}: Indices of subsystem sites to trace outstate::Vector{ET}: Full density matrix state in vectorized form
Returns
Matrix{Float64}: Reduced density matrix of the subsystem
For ladder systems where both subsystems have equal dimensions.
FibonacciChain.ladderbraidingsqmap — Function
ladderbraidingsqmap(model::AnyonModel{FibonacciAnyon}, state::Vector{ET}, idx::Int) where {ET}Apply braiding squared operation to density matrix state in vectorized form. Effectively, it describes the noise induced by inter-two layer chain braiding.
Arguments
model::AnyonModel{FibonacciAnyon}: Fibonacci anyon model containing system parameters (N, pbc)state::Vector{ET}: Density matrix state in vectorized formidx::Int: Site index for braiding operation
Returns
Vector{ET}: Transformed density matrix state after braiding
Operates on superposition states represented as vectorized density matrices.
Examples
julia> using FibonacciChain, LinearAlgebra, BitBasis
julia> N = 4;
julia> model = AnyonModel(FibonacciAnyon(), N; pbc=true);
julia> # Create PBC Fibonacci anyon basis
basis = anyon_basis(model);
julia> l = length(basis);
julia> ρ_vec = zeros(ComplexF64, l^2);
julia> # Initialize as identity matrix (vectorized)
for i in 1:l
ρ_vec[(i-1)*l + i] = 1.0/l;
end
julia> # Apply braiding at site 2
ρ_braided = ladderbraidingsqmap(model, ρ_vec, 2);
julia> norm(ρ_braided) > 0 # Should be non-zero
trueFibonacciChain.laddertranslationmap — Function
laddertranslationmap(model::AnyonModel, state::Vector{ET}) where {ET}Apply translation operator to vectorized density matrix state.
Arguments
model::AnyonModel: Anyon model containing system parameters (N, pbc, anyon_type)state::Vector{ET}: Density matrix state in vectorized form
Returns
Vector{ET}: Translated density matrix state
Translates both bra and ket parts of the density matrix consistently.
Measurement Functions
FibonacciChain.measure_basismap — Function
measure_basismap(model::AnyonModel, τ::Float64, state::T, i::Int, sign::Bool) where {T}Map single basis state under measurement operation at site i.
Arguments
model::AnyonModel: Anyon model containing system parameters (N, pbc, measure_operator)τ::Float64: Measurement strength parameterstate::T: Input basis statei::Int: Measurement site index (1 ≤ i ≤ N)sign::Bool: Measurement outcome (false for +, true for -)
Returns
Tuple: Either(basis1, basis2, coeff1, coeff2)for superposition output or(basis, coeff)for single output
Maps individual basis states according to measurement protocols and fusion rules.
Examples
julia> using FibonacciChain, BitBasis
julia> N = 4; T = BitStr{N, Int};
julia> model = AnyonModel(FibonacciAnyon(), N; pbc=true, measure_operator=:Antiferro);
julia> state = T(0b0100); # Single τ at site 2
julia> τ = 1.0; # Measurement strength
julia> result = measure_basismap(model, τ, state, 2, false);
julia> length(result) ∈ [2, 4] # Returns 2 or 4 elements depending on configuration
trueFibonacciChain.measuremap — Function
measuremap(model::AnyonModel, τ::Float64, state::Vector{ET}, idx::Int, sign::Bool)
measuremap(model::AnyonModel, ψ::MPS, sites, i::Int, τ::Float64, sign::Bool;
cutoff::Float64=1e-10, maxdim::Int=100)Apply measurement operator to quantum state and return post-measurement state.
Methods
Vector version (from Measurement.jl)
Apply measurement to a state vector.
Arguments
model::AnyonModel: Anyon model containing system parametersτ::Float64: Measurement strength parameterstate::Vector{ET}: Input quantum state vectoridx::Int: Measurement site indexsign::Bool: Measurement outcome (false for +, true for -)
Returns
Vector{ET}: Post-measurement quantum state (unnormalized)
MPS version (from MPSMeasurement.jl)
Apply measurement to an MPS state.
Arguments
model::AnyonModel: Anyon model containing system parametersψ::MPS: Input quantum statesites: ITensor site indicesi::Int: Measurement siteτ::Float64: Measurement strength parametersign::Bool: Measurement outcome (false for +, true for -)cutoff::Float64=1e-10: MPS truncation cutoffmaxdim::Int=100: Maximum bond dimension
Returns
MPS: Post-measurement quantum state (normalized)Float64: Measurement probability
FibonacciChain.measurement_enumeration — Function
measurement_enumeration(model::AnyonModel, τ::Float64, initial_state::Vector{ET}, measurement_sites::Vector{Int}) where {ET}Enumerate all trajectories of measurements on a given initial state.
Arguments
model::AnyonModel: Anyon model containing system parametersτ::Float64: Measurement strength parameterinitial_state::Vector{ET}: Initial quantum state vectormeasurement_sites::Vector{Int}: Sites to measure
Returns
Tuple{Vector{Vector{ET}}, Vector{Vector{Int64}}, Vector{Float64}}: (final_states, trajectories, probabilities)
Examples
julia> using FibonacciChain, LinearAlgebra
julia> model = AnyonModel(FibonacciAnyon(), 4; pbc=true, measure_operator=:Antiferro);
julia> basis = anyon_basis(model);
julia> state = normalize(ones(length(basis)));
julia> states, trajs, probs = measurement_enumeration(model, 1.0, state, [2]);
julia> length(states) == 2 # Two possible outcomes
trueFibonacciChain.measurement_tree_visualization — Function
measurement_tree_visualization(trajectories::Vector{Vector{Int64}}, probabilities::Vector{Float64})Visualize a measurement tree given trajectories and their probabilities.
Arguments
trajectories::Vector{Vector{Int64}}: Measurement outcome sequences at each branch.probabilities::Vector{Float64}: Corresponding probabilities for each trajectory.
Behavior
- Normalizes probabilities to sum to 1.
- Prints levels from root to leaves with indentation representing depth.
FibonacciChain.transfer_matrix — Function
transfer_matrix(model::AnyonModel, τ::Float64; sign::Bool=true)Construct the transfer matrix for two measurement layers.
Arguments
model::AnyonModel: Anyon model containing system parametersτ::Float64: Measurement strength parametersign::Bool=true: Measurement outcome sign
Returns
Matrix{Float64}: Transfer matrix between two measurement layers
FibonacciChain.boundary_evolution — Function
boundary_evolution(model::AnyonModel, state::Vector{T}, measure_config::MeasureConfig,
sample::Union{Nothing, BitVector}=nothing; layer_idx::Int=1)
boundary_evolution(model::AnyonModel, sites, state::MPS, sample::BitVector,
measure_config::MeasureConfig;
cutoff::Float64=1e-12, maxdim::Int=1000, layer_idx::Int=1)Evolve a quantum state under a single layer of boundary measurements with a given trajectory.
Methods
Vector version (from Measurement.jl)
Evolve a state vector under boundary measurements.
Arguments
model::AnyonModel: The anyon model specifying system parameters.state::Vector{T}: The initial state vector.measure_config::MeasureConfig: Configuration containing measurement strengthτand mode.sample::Union{Nothing, BitVector}=nothing: The measurement outcomes for the layer.layer_idx::Int=1: The layer index for measurement configuration.
Returns
Measurement_outcome_boundary: A struct containingstate,sample, andfree_energy.
MPS version (from MPSMeasurement.jl)
Evolve an MPS state under boundary measurements.
Arguments
model::AnyonModel: The anyon model specifying system parameters.sites: The ITensor site indices.state::MPS: The initial MPS state.sample::BitVector: The measurement outcomes for the layer.measure_config::MeasureConfig: Configuration containing measurement strengthτand mode (:sampleor:Born).cutoff::Float64=1e-12: MPS truncation cutoff.maxdim::Int=1000: Maximum bond dimension for MPS operations.layer_idx::Int=1: The layer index for measurement configuration.
Returns
Measurement_outcome_mps_boundary: A struct containing:state::MPS: The evolved MPS state.samples::BitVector: The measurement outcomes for the layer.free_energy::Float64: The free energy associated with the measurement layer.
FibonacciChain.bulk_evolution — Function
bulk_evolution(model::AnyonModel, state::Vector{ET}, measure_config::MeasureConfig,
samples::Union{Nothing,BitMatrix}=nothing)
bulk_evolution(model::AnyonModel, sites, state::MPS, measure_config::MeasureConfig,
samples::Union{Nothing,BitMatrix}=nothing;
cutoff::Float64=1e-10, maxdim::Int=100)Perform bulk measurement evolution from t₁ to t₂ on quantum state.
Methods
Vector version (from Measurement.jl)
Evolve a state vector under bulk measurements.
Arguments
model::AnyonModel: Anyon model containing system parametersstate::Vector{ET}: Initial quantum state vectormeasure_config::MeasureConfig: Configuration struct containing τ, t₁, t₂, mode, rng, etc.samples::Union{Nothing,BitMatrix}=nothing: Predefined measurement sequences for:samplemode
Returns
Measurement_outcome_bulk: A struct containing:states::Vector{Vector{ET}}: Intermediate states at each time stepsamples::BitMatrix: Measurement outcome sequencesfree_energys::Vector{Float64}: Free energy for each layer
MPS version (from MPSMeasurement.jl)
Evolve an MPS state under bulk measurements.
Arguments
model::AnyonModel: Anyon model containing system parameterssites: ITensor site indicesstate::MPS: Initial MPS quantum statemeasure_config::MeasureConfig: Configuration struct containing τ, t₁, t₂, mode, rng, etc.samples::Union{Nothing,BitMatrix}=nothing: Predefined measurement sequences for:samplemodecutoff::Float64=1e-10: MPS truncation cutoffmaxdim::Int=100: Maximum bond dimension
Returns
Measurement_outcome_bulk: A struct containing:states::Vector{MPS}: Intermediate states at each time stepsamples::BitMatrix: Measurement outcome sequencesfree_energy::Vector{Float64}: Free energy for each layer
Notes
- In
:Bornmode, samples are generated probabilistically via Born rule - In
:samplemode,samplesmust be provided as input - (2N+1) layers of measurements correspond to N time steps of evolution
MPS Functions
FibonacciChain.measurement_operator_mps — Function
measurement_operator_mps(model::AnyonModel, sites::Vector{<:Index}, i::Int, τ::Float64, sign::Bool)Create local measurement operator at site i as ITensor.
Arguments
model::AnyonModel: Anyon model containing system parameterssites::Vector{<:Index}: ITensor site indicesi::Int: Measurement siteτ::Float64: Measurement strength parametersign::Bool: Measurement outcome (false for +, true for -)
Returns
ITensor: Local measurement operator incorporating neighboring site correlations
FibonacciChain.initial_mps — Function
anyon_mps_gst(model::AnyonModel; sweep_times=5, maxdim=5, cutoff=1e-10, outputlevel=1)Find ground state of anyon chain Hamiltonian using DMRG.
Arguments
model::AnyonModel: Anyon model containing system parameters (N, pbc, anyon_type)sweep_times=5: Number of DMRG sweepsmaxdim=5: Maximum bond dimensioncutoff=1e-10: Truncation cutoffoutputlevel=1: Verbosity level
Returns
MPS: Ground state as Matrix Product StateFloat64: Ground state energy
Examples
julia> using FibonacciChain, ITensorMPS, ITensors
julia> model = AnyonModel(FibonacciAnyon(), 8; pbc=true);
julia> ψ_gs, E0 = anyon_mps_gst(model, maxdim=10, outputlevel=0);
julia> ψ_gs isa MPS
true
julia> E0 < 0 && imag(E0) ≈ 0
trueFibonacciChain.mps_measurement_enumeration — Function
mps_measurement_enumeration(model::AnyonModel, ψ::MPS, sites::Vector{<:Index}, measurement_sites::Vector{Int},
τ::Float64; cutoff::Float64=1e-10, maxdim::Int=100)Enumerate all possible measurement trajectories on MPS state.
Arguments
model::AnyonModel: Anyon model containing system parametersψ::MPS: Input quantum statesites: ITensor site indicesmeasurement_sites::Vector{Int}: Sites to measureτ::Float64: Measurement strength parametercutoff::Float64=1e-10: MPS truncation cutoffmaxdim::Int=100: Maximum bond dimension
Returns
Vector{MPS}: Final states for each trajectoryVector{Vector{Bool}}: Measurement trajectoriesVector{Float64}: Probabilities for each trajectory
Reference Qubit Functions
FibonacciChain.add_reference_qubits — Function
add_reference_qubits(model::AnyonModel, ψ::MPS, sites::Vector{<:Index}, site_idx::Int=1;
k_new::Int=1, verbose::Bool=false) -> Tuple{MPS, Vector{<:Index}}Add reference qubits to an MPS state at specified site using copy method for correlation measurements.
Principle
Reference qubits are ancillary qubits used to probe spatial and temporal correlations in monitored quantum systems. The key idea is to create a maximally entangled state between the reference qubit and a system qubit at a specific site, enabling measurement of correlations through the reference qubit's reduced density matrix.
Copy Method
Creates a controlled-copy (CNOT-like) entanglement:
- If the system qubit at
site_idxis in state |0⟩, the reference qubit becomes |0⟩ - If the system qubit at
site_idxis in state |1⟩, the reference qubit becomes |1⟩
Mathematically, for a state |ψ⟩ = α|0⟩ + β|1⟩ at site_idx:
|ψ⟩ → α|00⟩ + β|11⟩ (reference ⊗ system)This preserves the original state's structure while creating classical correlation.
For reset method (measure system qubit first, then create Bell pair), use add_reference_qubits_reset instead.
Arguments
model::AnyonModel: Anyon model containing system parametersψ::MPS: Input MPS quantum statesites: ITensor site indicessite_idx::Int=1: Site index to entangle with reference qubit (1 ≤ site_idx ≤ N)k_new::Int=1: Number of reference qubits to add (must be 0 or 1)verbose::Bool=false: Enable verbose output for debugging
Returns
MPS: New MPS with reference qubit added at the leftmost positionVector{Index}: Updated site indices including the new reference qubit
Notes
- Each system qubit can only be entangled with one reference qubit
- To add multiple reference qubits at different sites, call this function multiple times
- The reference qubit is inserted at position 1 (leftmost) of the MPS
- Original MPS is not modified (returns a new copy)
Example
model = AnyonModel(FibonacciAnyon(), 6)
ψ, sites = anyon_mps_gst(model)
ψ_ref, sites_ref = add_reference_qubits(model, ψ, sites, 3) # Add ref qubit at site 3See also: add_reference_qubits_reset
add_reference_qubits(model::AnyonModel, state::Vector{ET}, site_idx::Int64;
k_new::Int=1, verbose::Bool=false) -> Vector{ET}Add reference qubits to a quantum state at specified site using copy method for correlation measurements.
Principle
Reference qubits are ancillary qubits used to probe spatial and temporal correlations in monitored quantum systems. The key idea is to create a maximally entangled state between the reference qubit and a system qubit at a specific site, enabling measurement of correlations through the reference qubit's reduced density matrix.
Copy Method
Creates a controlled-copy (CNOT-like) entanglement:
- If the system qubit at
site_idxis in state |0⟩, the reference qubit becomes |0⟩ - If the system qubit at
site_idxis in state |1⟩, the reference qubit becomes |1⟩
Mathematically, for a state |ψ⟩ = α|0⟩ + β|1⟩ at site_idx:
|ψ⟩ → α|00⟩ + β|11⟩ (system ⊗ reference)
|ψ⟩ → α|A0⟩ + β|B1⟩ (system ⊗ reference), where A, B are orthogonal statesThis preserves the original state's structure while creating classical correlation.
For reset method (measure system qubit first, then create Bell pair), use add_reference_qubits_reset instead.
Basis Structure
The extended basis combines reference qubits (prefix) with the Fibonacci/Ising constraint basis (suffix):
|extended_basis⟩ = |ref₁ ref₂ ... refₖ⟩ ⊗ |system basis⟩Reference qubits are stored as the leftmost (highest) bits in the binary representation.
Arguments
model::AnyonModel: Anyon model containing system size N, boundary conditions, etc.state::Vector{ET}: Quantum state vector (may already contain k_old reference qubits)site_idx::Int64: Site index for reference qubit insertion (1 ≤ site_idx ≤ N)k_new::Int=1: Number of new reference qubits to add (must be 0 or 1)verbose::Bool=false: Enable verbose output for debugging
Returns
Vector{ET}: New normalized state with reference qubit added
Notes
- Each system qubit can only be entangled with one reference qubit
- To add multiple reference qubits at different sites, call this function multiple times
- The number of existing reference qubits (k_old) is automatically deduced from state length
- State dimension grows as: newdim = 2^(kold + knew) × len(Fibonaccibasis)
Example
# Add first reference qubit at site 1
state1 = add_reference_qubits(model, initial_state, 1)
# Add second reference qubit at site N÷2
state2 = add_reference_qubits(model, state1, N÷2)
# Now state2 has 2 reference qubits for 2-point correlation measurementSee also: add_reference_qubits_reset, reference_rdm, reference_evolution, build_extended_basis
FibonacciChain.add_reference_qubits_reset — Function
add_reference_qubits_reset(model::AnyonModel, ψ::MPS, sites::Vector{<:Index}, site_idx::Int=1;
k_new::Int=1, verbose::Bool=false) -> Tuple{Float64, Float64, MPS, MPS, Vector{Index}}Add reference qubits using reset method: first measures/resets the system qubit, then creates a Bell pair.
Principle
- Project the system qubit onto |0⟩ or |1⟩ basis (strong measurement)
- Create a Bell pair between the reset qubit and reference qubit
- Returns both branches with their probabilities
This method is useful when you want to prepare a known initial state before adding the reference qubit.
Arguments
model::AnyonModel: Anyon model containing system parametersψ::MPS: Input MPS quantum statesites: ITensor site indicessite_idx::Int=1: Site index for reference qubit insertion (1 ≤ site_idx ≤ N)k_new::Int=1: Number of new reference qubits to add (must be 0 or 1)verbose::Bool=false: Enable verbose output for debugging
Returns
Float64: Probability of measuring |0⟩ at site_idxFloat64: Probability of measuring |1⟩ at site_idxMPS: State after measuring |0⟩ and adding Bell pairMPS: State after measuring |1⟩ and adding Bell pairVector{Index}: Updated site indices including the new reference qubit
Notes
- Unlike
add_reference_qubits, this method returns both measurement branches - The probabilities satisfy prob₀ + prob₁ = 1
- Each returned MPS is normalized
Example
model = AnyonModel(FibonacciAnyon(), 6)
ψ, sites = anyon_mps_gst(model)
prob0, prob1, ψ0, ψ1, sites_ref = add_reference_qubits_reset(model, ψ, sites, 3)
# Use ψ0 with probability prob0, or ψ1 with probability prob1See also: add_reference_qubits
add_reference_qubits_reset(model::AnyonModel, state::Vector{ET}, site_idx::Int64;
k_new::Int=1, verbose::Bool=false) -> Tuple{Float64, Float64, Vector{ET}, Vector{ET}}Add reference qubits using reset method: first measures/resets the system qubit, then creates a Bell pair.
Principle
- Project the system qubit onto |0⟩ or |1⟩ basis (strong measurement)
- Create a Bell pair between the reset qubit and reference qubit
- Returns both branches with their probabilities
This method is useful when you want to prepare a known initial state before adding the reference qubit.
Arguments
model::AnyonModel: Anyon model containing system size N, boundary conditions, etc.state::Vector{ET}: Quantum state vector (may already contain k_old reference qubits)site_idx::Int64: Site index for reference qubit insertion (1 ≤ site_idx ≤ N)k_new::Int=1: Number of new reference qubits to add (must be 0 or 1)verbose::Bool=false: Enable verbose output for debugging
Returns
Tuple{Float64, Float64, Vector{ET}, Vector{ET}}: (prob₀, prob₁, stateafter0, stateafter1) for both measurement branches
See also: add_reference_qubits, reference_rdm, reference_evolution
FibonacciChain.reference_measuremap — Function
reference_measuremap(model::AnyonModel, τ::Float64, state::Vector{ET}, idx::Int, sign::Bool;
extended_basis::Vector{newT}, k_old::Int64=1)
reference_measuremap(model::AnyonModel, ::Type{T}, ::Type{pretype}, τ::Float64, state::Vector{ET},
idx::Int, sign::Bool; extended_basis::Vector{newT})Apply a measurement operator to a quantum state with reference qubits.
This function extends measuremap to handle states that include reference qubits, preserving the reference qubit subspace while applying measurements to the system.
Arguments
model::AnyonModel: Anyon model containing system parametersτ::Float64: Measurement strength parameterstate::Vector{ET}: Quantum state vector with reference qubitsidx::Int: Site index for measurementsign::Bool: Measurement outcome (false=0, true=1)extended_basis::Vector{newT}: Extended basis including reference qubitsk_old::Int64=1: Number of reference qubits in the state
Returns
Vector{ET}: The post-measurement state (unnormalized)
Notes
The extended basis must be constructed via build_extended_basis before calling this function.
FibonacciChain.reference_boundary_evolution — Function
reference_boundary_evolution(model::AnyonModel, state::Vector{T}, measure_config::MeasureConfig,
sample::Union{Nothing, BitVector}=nothing; layer_idx::Int64=1)Perform a single measurement layer evolution on a state with reference qubits.
Arguments
model::AnyonModel: Anyon model containing system parametersstate::Vector{T}: Quantum state vector with reference qubitsmeasure_config::MeasureConfig: Configuration containing τ, mode, rng, etc.sample::Union{Nothing, BitVector}=nothing: Predefined measurement outcomes for:samplemodelayer_idx::Int64=1: Layer index (1-based) to determine measurement pattern
Returns
Measurement_outcome_boundary: A struct containing:state::Vector{T}: Post-measurement statesample::BitVector: Measurement outcomes for this layerfree_energy::Float64: Free energy contribution from this layer
Notes
- In
:Bornmode, samples are generated probabilistically via Born rule - In
:samplemode,samplemust be provided as input
FibonacciChain.reference_bulk_evolution — Function
reference_bulk_evolution(model::AnyonModel, state::Vector{T}, measure_config::MeasureConfig,
samples::Union{Nothing,BitMatrix}=nothing) where{T}Perform bulk measurement evolution from t₁ to t₂ on quantum state with reference qubits.
Arguments
model::AnyonModel: Anyon model containing system parametersstate::Vector{T}: Initial quantum state vector with reference qubitsmeasure_config::MeasureConfig: Configuration struct containing τ, t₁, t₂, mode, rng, etc.samples::Union{Nothing,BitMatrix}=nothing: Predefined measurement sequences for:samplemode
Returns
Measurement_outcome_bulk: A struct containing:states::Vector{Vector{T}}: Intermediate states at each time stepsamples::BitMatrix: Measurement outcome sequencesfree_energys::Vector{Float64}: Free energy for each layer
Notes
- In
:Bornmode, samples are generated probabilistically via Born rule - In
:samplemode,samplesmust be provided as input - The state should already contain reference qubits added via
add_reference_qubits
FibonacciChain.reference_rdm — Function
reference_rdm(model::AnyonModel, subsystems::Vector{Int64}, state::Vector{ET}; traceref::Bool=true)Compute the reduced density matrix for a state with reference qubits.
Arguments
model::AnyonModel: Anyon model containing system parameterssubsystems::Vector{Int64}: Indices of system sites to keep (not trace out)state::Vector{ET}: Quantum state vector with reference qubitstraceref::Bool=true: Iftrue, trace out reference qubits and keep system subsystems; iffalse, trace out system and keep reference qubits
Returns
Matrix{ET}: The reduced density matrix
Notes
- Subsystem indices count from the right of the binary string representation
- The number of reference qubits
k_oldis automatically deduced from the state length
FibonacciChain.reference_evolution — Function
reference_evolution(N::Int, τ::Float64, forward::Vector{ET}, sample::Matrix{Bool},
x₂::Int, t₁, t₂; x₁::Int=1, rng=MersenneTwister(), pbc=true,
anyon_type::Symbol=:Fibo, verbose=false, mode::Symbol=:sample)
reference_evolution(model::AnyonModel, sites, forward::Vector{MPS}, sample::Matrix{Bool},
x₂::Int, t₁, t₂, measure_config::MeasureConfig; x₁::Int=1, verbose=false)Compute temporal/spatial correlation between two spacetime points using cached forward evolution.
This function avoids recomputing the forward evolution up to time t₁ by using the cached forward states.
Methods
Vector version (from ReferenceProbe.jl)
Compute correlation using state vectors.
Arguments
N::Int: System sizeτ::Float64: Measurement strength parameterforward::Vector{ET}: Cached forward state evolution trajectorysample::Matrix{Bool}: Measurement sample configurationx₂::Int: Site index for second reference qubit insertiont₁::Int: First time slice indext₂::Int: Second time slice index (must be >= t₁)x₁::Int=1: Spatial site index for first reference qubitrng::MersenneTwister=MersenneTwister(): Random number generatorpbc::Bool=true: Periodic boundary conditionsanyon_type::Symbol=:Fibo: Anyon type (:Fiboor:Ising)verbose::Bool=false: Enable verbose outputmode::Symbol=:sample: Evolution mode (:sampleor:Born)
Returns
- Reference qubit added states between specified spacetime slices.
MPS version (from MPSMeasurement.jl)
Compute correlation using MPS states.
Arguments
model::AnyonModel: Anyon model containing system parameterssites: ITensor site indicesforward::Vector{MPS}: Cached forward state evolution trajectory from a previousbulk_evolutionrun.sample::Matrix{Bool}: The full measurement sample configuration for the entire evolution.x₂::Int: The site index for the second reference qubit insertion.t₁::Int: The first time slice index for correlation.t₂::Int: Second time slice index (must be >= t₁)measure_config::MeasureConfig: Configuration containing τ, rng, mode, etc.x₁::Int=1: Site index for first reference qubitverbose::Bool=false: Enable verbose output for debugging.
Returns
Vector{MPS}: A vector of MPS states for the specified spacetime slices.BitMatrix: The measurement samples used for the evolution.Vector{Float64}: The free energy calculated for each layer.
reference_evolution(model::AnyonModel, forward::Vector{ET}, measure_config::MeasureConfig,
sample::Union{BitMatrix, Nothing}=nothing) -> Measurement_outcome_bulkCompute temporal/spatial correlation between two spacetime points using cached forward evolution.
This function avoids recomputing the forward evolution up to time t₁ by using the cached forward states, then adds reference qubits and continues evolution to measure correlations.
Spacetime Diagram
| ---->|
forward
Ref1 x₁ = 1
|
|
Ref2 --> Ref3 x₂ = N/2
| ----> |______| -----> |
1 t₁ t₂=t₁+δt Δt
| --------------------->|
sampleArguments
model::AnyonModel: Anyon model containing system parameters (N, pbc, anyon_type, etc.)forward::Vector{ET}: Cached forward state evolution trajectory from a previousbulk_evolutionrunmeasure_config::MeasureConfig: Configuration struct containing:τ::Float64: Measurement strength parametert₁::Int: First time slice index for correlationt₂::Int: Second time slice index (must be >= t₁)x₁::Int: Site index for first reference qubitx₂::Int: Site index for second reference qubit (must be >= x₁)mode::Symbol: Evolution mode (:sampleor:Born)rng::MersenneTwister: Random number generator (for:Bornmode)verbose::Bool: Enable verbose output
sample::Union{BitMatrix, Nothing}=nothing: Full measurement sample configuration for the entire evolution
Returns
Measurement_outcome_bulk: A struct containing:states::Vector{ET}: State evolution trajectory with reference qubitssamples::BitMatrix: Measurement samples used for the evolutionfree_energys::Vector{Float64}: Free energy calculated for each layer
Correlation Types
Based on δt = t₂ - t₁ and δx = |x₂ - x₁|:
δt > 0 && δx > 0: 3-point correlation (both spatial and temporal), uses 3 reference qubitsδt == 0: Pure 2-point spatial correlation at fixed time, uses 2 reference qubitsδx == 0: Pure 2-point temporal correlation at fixed site, uses 2 reference qubits
Notes
- The
forwardstates should come from a previousbulk_evolutioncall - Reference qubits are added via
add_reference_qubitsat the specified spacetime points - In
:samplemode,samplemust be provided; in:Bornmode, samples are generated via Born rule
See also: add_reference_qubits, reference_bulk_evolution, ref_correlation