Documentation

Basics

PauliStrings.OperatorTSMethod
OperatorTS{Ls}(o::Operator)
OperatorTS{Ls,Ps}(o::Operator)

Construct an $n$-dimensional translation symmetric operator from o where Ls is a tuple of integers (L1, L2, ...) and Ps is an optional tuple of Bool values indicating which dimensions are periodic (default: all true). The resulting operator is equivalent to

\[O_\mathrm{TS} = \sum_T T^\dag O T\]

where $T$ are all translations on the L1×L2×… hypercube. So if you feed it an operator that is already a sum, you should afterwards normalize it by the number of sites.

To get a dense operator from this lazy sum representation, see resum. To get a single term, see representative.

source
PauliStrings.qubitlengthFunction
qubitlength(x::AbstractOperator)
qubitlength(::Type{<:AbstractOperator})

Returns the number of qubits the operator acts on.

source

Truncation

Base.truncateMethod
truncate(o::Operator, max_lenght::Int; keepnorm::Bool = false)

Remove all terms of length > maxlenght. Keep all terms of length <= maxlenght. i.e remove all M-local terms with M>max_lenght

Example

A = Operator(4)
A += "X",1,"X",2
A += "Z",1,"Z",2,"Z",4
julia> A
(1.0 + 0.0im) ZZ1Z
(1.0 + 0.0im) XX11

julia> ps.truncate(A,2)
(1.0 + 0.0im) XX11
source
PauliStrings.k_local_partFunction
k_local_part(o::Operator, k::Int; atmost=false)

Return the k-local part of o. I.e all the strings of lenght k. Set atmost=true to return the 'at most' k-local part, i.e all the strings of length <= k.

Example

A = Operator(4)
A += "X",1,"X",2
A += "Z",1,"Z",2,"Z",4
A += "1X11"
julia> A
(1.0 + 0.0im) ZZ1Z
(1.0 + 0.0im) 1X11
(1.0 + 0.0im) XX11

julia> k_local_part(A,2)
(1.0 + 0.0im) XX11
source
PauliStrings.trimFunction
trim(o::Operator, max_strings::Int; keepnorm::Bool = false, keep::Operator=Operator(N))

Keep the first max_strings terms with largest coeficients.

keepnorm is set to true to keep the norm of o.

keep is an operator that specify a set of strings that cannot be removed

Example

A = Operator(4)
A += 1,"XXXX"
A += 2,"XX11"
A += 3,"XX1X"
A += 4,"ZZXX"
B = Operator(4)
B += 1,"XX11"
B += 1,"XX1X"
julia> trim(A,2)
(4.0 + 0.0im) ZZXX
(3.0 + 0.0im) XX1X

julia> trim(A,2;keep=B)
(4.0 + 0.0im) ZZXX
(3.0 + 0.0im) XX1X
(2.0 + 0.0im) XX11
source
PauliStrings.pruneFunction
 prune(o::Operator, alpha::Real; keepnorm::Bool = false)

Keep terms with probability 1-exp(-alpha*abs(c)) where c is the weight of the term

source
PauliStrings.cutoffFunction
cutoff(o::Operator, epsilon::Real; keepnorm::Bool = false)

Remove all terms with weight < epsilon

Example

A = Operator(4)
A += 1,"XXXX"
A += 2,"XX11"
A += 3,"XX1X"
A += 4,"ZZXX"
julia> cutoff(A, 2.5)
(3.0 + 0.0im) XX1X
(4.0 + 0.0im) ZZXX
source

Noise

PauliStrings.add_noiseFunction
add_noise(o::Operator, g::Real)

Add depolarizing noise that make the long string decays. g is the noise amplitude. Each string is multiplied by exp(-g * w), where w is the number of non-identity Pauli operators in the string. This is equivalent to the compostion of $N$ single qubit depolarizing channels with Kraus operators $\sqrt{1-\frac{3p}{4}} I_i, \, \sqrt{\frac{p}{4}} X_i, \, \sqrt{\frac{p}{4}} Y_i, \, \sqrt{\frac{p}{4}} Z_i$ and $p=1-e^{-g}$.

Example

A = add_noise(A, 0.1)

Reference

https://arxiv.org/pdf/2407.12768

source
add_noise(o::Operator, g::AbstractVector{<:Real})

Add local depolarizing noise.

If $g_j$ is the noise amplitude for site $j$, then each string will be multiplied by $e^{-\sum_j g_j}$, where the sum runs over the sites with non-unit Pauli operators.

source
PauliStrings.add_dephasing_noiseFunction
add_dephasing_noise(o::AbstractOperator, g::Real; basis::Symbol=:Z)

Add dephasing noise.

If $g$ is the noise amplitude, then each string will decay by a factor of $e^{-gw}$, where $w$ is the count of Pauli operators in the string that are either $X$ or $Y$.

Reference

https://arxiv.org/abs/2306.05804

source
add_dephasing_noise(o::AbstractOperator, g::AbstractVector{<:Real})

Add local dephasing noise.

If $g_j$ is the noise amplitude of site $j$, then each string will be multiplied by $e^{-\sum_j g_j}$, where the sum runs over the sites with Pauli operators that are either $X$ or $Y$.

source

Time evolution

PauliStrings.evolveFunction
evolve(H::AbstractOperator, O::AbstractOperator, tspan::AbstractVector; kwargs...) -> EvolutionResult
evolve(H::AbstractOperator, O::AbstractOperator, dt::Real, nsteps::Integer; kwargs...) -> EvolutionResult

Evolve the operator O under the Hamiltonian H in the Heisenberg picture and return an EvolutionResult. Values are saved at the times in tspan; tspan[1] is saved before any integration.

The integrator takes one internal step per save-interval. To use a finer internal step than the spacing at which results are saved, pass a finer tspan.

Keyword arguments

  • method::AbstractEvolutionMethod = RK4(). One of Trotter, RK4, DOPRI5, Exact.
  • truncation: function O -> O applied after every internal step. Default identity.
  • dissipation: function (O, dt) -> O applied after every internal step. The dt argument reflects that dissipation is part of the dynamics (a Lindblad semigroup factor) and depends on the step size. Default (O, dt) -> O.
  • fout: function fout(O) called at every save time; its return value is collected into history. Default (O) -> nothing. Pass fout = copy to save the full operator trajectory. Pass fout = nothing to skip history collection entirely. final is populated independently of fout.
  • hbar::Real = 1: Planck constant, forwarded to the underlying step routine.

Examples

Full-trajectory (fout = copy):

result = evolve(H, O, 0.0:0.05:1.0; fout = copy)
result.history          # Vector of operators saved at each time
result.final            # final operator

Scalar observable trajectory:

result = evolve(H, O, 0.0:0.05:1.0; fout = O -> trace(O * A))
result.history          # Vector{ComplexF64} of ⟨A(t)⟩
result.final            # final operator, free for continuation

Skip history (default), keep only the final state:

result = evolve(H, O, 0.0:0.05:1.0; fout = nothing)
result.history === nothing
result.final            # final operator

Convenience (dt, nsteps) form:

evolve(H, O, 0.01, 100)        # 100 steps of dt = 0.01, 101 save points

Dissipative run:

evolve(H, O, 0.0:0.05:1.0; dissipation = (O, dt) -> add_noise(O, 0.1 * dt))

Choose an integrator:

evolve(H, O, 0.0:0.05:1.0; method = RK4())
evolve(H, O, 0.0:0.05:1.0; method = DOPRI5())
source
PauliStrings.EvolutionResultType
EvolutionResult{T, H, O}

Result of evolve.

Fields

  • t::Vector{T}: the save times (collect(tspan)).
  • history::H: a Vector whose i-th element is fout(O(t[i]), t[i]), or nothing if fout === nothing was passed.
  • final::O: the final operator, always populated regardless of fout.
source
PauliStrings.AbstractEvolutionMethodType
AbstractEvolutionMethod

Supertype for time-evolution algorithm specifications passed to evolve. Each concrete subtype carries only that algorithm's own options. New methods are new subtypes; the signature of evolve does not change.

source
PauliStrings.RK4Type
RK4()

Classical fixed-step 4th-order Runge–Kutta. Takes one internal step per save-interval; subdivide by supplying a finer tspan.

source
PauliStrings.DOPRI5Type
DOPRI5()

Dormand–Prince 5(4) embedded Runge–Kutta. In this version it is used as a fixed-step method: one internal step per save-interval.

source
PauliStrings.TrotterType
Trotter(; order=2, gates=nothing)

Fixed-step product-formula integrator. order is 1 (Lie) or 2 (Strang). gates is an optional precomputed gate list (from trotterize) matching the save-interval dt; if nothing, gates are built internally from H and cached when tspan is uniformly spaced.

Implemented for H::Operator/O::Operator and for H::OperatorTS/O::OperatorTS. In the translation-symmetric case the gates are built from resum(H) and applied to the representative term, so a precomputed gates list must come from trotterize(resum(H), dt).

source
PauliStrings.ExactType
Exact()

Reference integrator. Builds dense matrices Hm = Matrix(H) and Om = Matrix(O) once, then evolves with Om <- U * Om * U' where U = exp(i Hm dt / hbar), exponentiating per save-interval (cached when tspan is uniform).

history elements and final are dense matrices, not Operators. Intended for testing/benchmarking small systems against the truncated-Pauli methods; cost and memory are exponential in the qubit count.

source
PauliStrings.TrotterGateType
TrotterGate{P,T}

One factor in a product-form Trotter step: on each Pauli string, the Heisenberg update exp(im * theta * G/2) * P * exp(-im * theta * G/2) with generator G, implemented via pauli_rotation. The theta field follows that convention.

source
PauliStrings.trotterizeFunction
trotterize(H::Operator, dt::Real; order=2, heisenberg=true, hbar=1)

Build a first-order (order=1, Lie) or second-order (order=2, Strang) Trotter list that approximates exp(im * H * dt / hbar) (Heisenberg) or the conjugate sequence (Schrödinger / density matrix). Each gate uses pauli_rotation with the returned theta field.

For H::Operator{<:PauliStringTS}, see the specialized trotterize that calls resum first.

source
PauliStrings.trotter_step!Function

trotterstep!(O::Operator, gates::AbstractVector{<:TrotterGate}; truncation::Function=identity, truncateevery::Int=1)

Apply one Trotter step in place. Gates must be listed in matrix-multiply order U = V1 * V2 * ... * Vn; conjugation O -> U * O * U' applies factors Vn, ..., V1 successively (reverse of the list). Each Pauli string uses the same coefficient convention as Matrix(O) (weights include division by im to the number of Y factors).

source
PauliStrings.pauli_rotationFunction
pauli_rotation(G::PauliString, P::PauliString, theta::Real)

Heisenberg-picture action of a Pauli rotation generated by G on a single Pauli string P, corresponding to

RG(theta)[P] = exp(1im * theta * G / 2) * P * exp(-1im * theta * G / 2)

as in Eqs. (13)–(14) of https://arxiv.org/pdf/2505.21606. If P commutes with G, this returns P as an Operator. Otherwise it returns the 2-branching combination cos(theta) P + sin(theta) P′ as an Operator, where P′ = (im/2) * [G, P] is constructed from the commutator.

source
PauliStrings.rk4Function
rk4(H::AbstractOperator, O::AbstractOperator, dt::Real; hbar::Real=1, heisenberg=true, truncation::Function=copy, f=f_unitary)

Single step of Runge–Kutta-4 with time independant Hamiltonian. Returns O(t+dt). truncation : function that takes an operator and returns a truncated version of it. By default it is identity (no truncation).

source
rk4(H::Function, O::AbstractOperator, dt::Real, t::Real; hbar::Real=1, heisenberg=true, truncation::Function=identity, f=f_unitary)

Single step of Runge–Kutta-4 with time dependant Hamiltonian. Returns O(t+dt). H is a function that takes time t and returns the Hamiltonian at that time. truncation : function that takes an operator and returns a truncated version of it. By default it is identity (no truncation).

source
PauliStrings.rk4_lindbladFunction
rk4_lindblad(H::AbstractOperator, O::AbstractOperator, dt::Real, L; hbar::Real=1, heisenberg=true, truncation::Function=copy, gamma=[]))

Single step of Runge–Kutta-4 for solving the Lindblad equation

$\dot{O}=i[H,O]+\sum_i \gamma_i \left(L_i^\dagger O L_i -\frac{1}{2} \{ L_i^\dagger L_i, O\} \right)$

Returns O(t+dt). L is a list of jump operators.

source

Other algorithms

PauliStrings.lanczosFunction
lanczos(H::Operator, O::Operator, steps::Int, nterms::Int; keepnorm=true, maxlength=1000, returnOn=false)

Compute the first steps lanczos coeficients for Hamiltonian H and initial operator O

At every step, the operator is trimed with trim and only nterms are kept.

Using maxlength speeds up the commutator by only keeping terms of length <= maxlength

Set returnOn=true to save the On's at each step. Then the function returns a pair of lists (bn, On). The first operators of the list On is O

source
PauliStrings.liomsFunction
lioms(H::AbstractOperator, support::Vector{T}; threshold::Real=1e-14, f::Function=(H,O)->im * commutator(H,O)) where {T<:AbstractOperator}

Algorithm constructing all Local Integrals of Motion (LIOMs) for a Hamiltonian H, supported on the operator basis from support. Follows definitions in https://arxiv.org/abs/2505.05882.

Arguments

  • H::T: The Hamiltonian operator
  • support::Vector{T}: Vector of basis operators
  • threshold::Real=1e-14: Threshold for eigenvalues above which eigenmodes are discarded. By default only exact LIOMs are returned.
  • f::Function=f: Function to check for LIOMs, by default the commutator f(H,O) = im*[H,O]

Returns

  • evals::Vector{Float64}: Eigenvalues below threshold
  • evecs::Matrix{Float64}: Eigenvectors corresponding to eigenvalues
  • ops::Vector{AbstractOperator}: LIOM operators
source
lioms(H::T, k::Int; threshold::Real=1e-14, f::Function=(H,O)->im * commutator(H,O)) where {T<:AbstractOperator}

Algorithm constructing all Local Integrals of Motion (LIOMs) for a Hamiltonian H, supported on the most general Pauli string basis on k sites. Follows definitions in https://arxiv.org/abs/2505.05882.

Arguments

  • H::T: The Hamiltonian operator
  • k::Int: Maximum support size for basis operators
  • threshold::Real=1e-14: Threshold for eigenvalues above which eigenmodes are discarded. By default only exact LIOMs are returned.
  • f::Function=f: Function to check for LIOMs, by default the commutator f(H,O) = im*[H,O]

Returns

  • evals::Vector{Float64}: Eigenvalues below threshold
  • evecs::Matrix{Float64}: Eigenvectors corresponding to eigenvalues
  • ops::Vector{T}: LIOM operators
source

Operations

Base.:+Method
Base.:+(o1::O, o2::O) where {O<:AbstractOperator}
Base.:+(o::AbstractOperator, a::Number)
Base.:+(a::Number, o::AbstractOperator)

Add two operators together or add a number to an operator

Example

A = Operator(4)
A += "XYZ1"
A += 1, "Y", 4
B = Operator(4)
B += 2, "Y", 2, "Y", 4
B += 1, "Z", 3
julia> A
(1.0 - 0.0im) 111Y
(1.0 - 0.0im) XYZ1

julia> B
(1.0 + 0.0im) 11Z1
(2.0 - 0.0im) 1Y1Y

julia> A+B
(1.0 + 0.0im) 11Z1
(2.0 - 0.0im) 1Y1Y
(1.0 - 0.0im) 111Y
(1.0 - 0.0im) XYZ1

julia> A+5
(1.0 - 0.0im) 111Y
(1.0 - 0.0im) XYZ1
(5.0 + 0.0im) 1111
source
Base.:-Method
Base.:-(o1::O, o2::O) where {O<:AbstractOperator}
Base.:-(o::AbstractOperator)
Base.:-(o::AbstractOperator, a::Number)
Base.:-(a::Number, o::AbstractOperator)
Base.:-(o1::Operator, o2::Operator)

Subtraction between operators and numbers

source
Base.:*Method
Base.:*(o1::Operator, o2::Operator; kwargs...)
Base.:*(o::Operator, a::Number)
Base.:*(a::Number, o::AbstractOperator)

Multiply two operators together or an operator with a number

Example

A = Operator(4)
A += "XYZ1"
A += 1, "Y", 4
B = Operator(4)
B += 2, "Y", 2, "Y", 4
B += 1, "Z", 3
julia> A
(1.0 - 0.0im) 111Y
(1.0 - 0.0im) XYZ1


julia> B
(1.0 + 0.0im) 11Z1
(2.0 - 0.0im) 1Y1Y

julia> A*B
(2.0 - 0.0im) X1ZY
(1.0 - 0.0im) 11ZY
(2.0 - 0.0im) 1Y11
(1.0 - 0.0im) XY11

julia> A*5
(5.0 - 0.0im) 111Y
(5.0 - 0.0im) XYZ1
source
Base.:^Method
Base.:^(o::Operator, k::Int)

kth power of o.

source
PauliStrings.commutatorFunction
commutator(o1::Operator, o2::Operator; kwargs...)

Commutator of two operators. This is faster than doing o1*o2 - o2*o1.

Example

julia> A = Operator(4)
julia> A += "X111"
julia> B = Operator(4)
julia> B += "Z111"
julia> B += "XYZ1"
julia> commutator(A,B)
(0.0 - 2.0im) Y111
source
PauliStrings.anticommutatorFunction
anticommutator(o1::Operator, o2::Operator; kwargs...)

Commutator of two operators. This is faster than doing o1*o2 + o2*o1.

source
Base.:/Method
Base.:/(o::AbstractOperator, a::Number)

Divide an operator by a number

source
PauliStrings.traceMethod
trace(o::Operator; normalize=false)
trace(o::OperatorTS)

Trace of an operator. If normalize is true, return the trace divided by 2^N.

Example

julia> A = Operator(4)
julia> A += 2,"1111"
julia> A += 3,"XYZ1"
julia> trace(A)
32.0 + 0.0im
source
LinearAlgebra.diagMethod
LinearAlgebra.diag(o::AbstractOperator)

Diagonal of an operator. Keep the strings that only contain 1's or Z's. Return another operator.

Example

julia> A = Operator(4)
julia> A += 2,"1111"
julia> A += 3,"XYZ1"
julia> A += 3,"Z11Z"
julia> diag(A)
(2.0 + 0.0im) 1111
(3.0 + 0.0im) Z11Z
source
LinearAlgebra.normMethod
LinearAlgebra.norm(o::AbstractOperator; normalize=false)

Frobenius norm, equivalent to sqrt(trace(o' * o)). If normalize is true, divide by sqrt(2^N).

Example

julia> A = Operator(4)
julia> A += 2,"X",2
julia> A += 1,"Z",1,"Z",3
julia> norm(A)
8.94427190999916
source
Base.adjointMethod
Base.adjoint(o::AbstractOperator)

Conjugate transpose. ' also works.

Example

A = Operator(3)
A += 1im,"X",2
A += 1,"Z",1,"Z",3
julia> A

(1.0 + 0.0im) Z1Z
(0.0 + 1.0im) 1X1


julia> adjoint(A)
(1.0 - 0.0im) Z1Z
(0.0 - 1.0im) 1X1

julia> A'
(1.0 - 0.0im) Z1Z
(0.0 - 1.0im) 1X1
source
PauliStrings.ptraceMethod
ptrace(o::AbstractOperator, keep::Vector{Int})

Partial trace.

keep is list of qubits indices to keep starting at 1 note that this still returns an operator of size N and doesnt permute the qubits this only gets rid of Pauli strings that have no support on keep and add their coeficient*2^N to the identity string

Example

A = Operator(5)
A += "XY1XZ"
A += "XY11Z"
julia> ptrace(A, [3,4])
(1.0 - 0.0im) XY1XZ
(8.0 - 0.0im) 11111

julia> ptrace(A, [1,5])
(1.0 - 0.0im) XY1XZ
(1.0 - 0.0im) XY11Z
source

Power and moments

Base.:^Method
Base.:^(o::Operator, k::Int)

kth power of o.

source
PauliStrings.trace_productFunction
trace_product(o1::Operator, o2::Operator; scale=0)
trace_product(o1::OperatorTS, o2::OperatorTS; scale=0)

Efficiently compute trace(o1*o2). This is much faster than doing trace(o1*o2). If scale is not 0, then the result is normalized such that trace(identity)=scale.

source
trace_product(A::Operator, k::Int, B::Operator, l::Int; scale=0)

Efficiently compute trace(A^k*B^l). This is much faster than doing trace(A^k*B^l).

If scale is not 0, then the result is normalized such that trace(identity)=scale.

source
trace_product(A::AbstractOperator; scale=0)

Compute trace(A*A). This is much faster than doing trace(A*A).

If scale is not 0, then the result is normalized such that trace(identity)=scale.

source
trace_product(A::Operator{<:PauliStringTS}; scale=0)

Compute trace(A*A). This is much faster than doing trace(A*A).

If scale is not 0, then the result is normalized such that trace(identity)=scale.

source
trace_product(A::AbstractOperator, k::Int; scale=0)

Efficiently compute trace(A^k). This is much faster than doing trace(A^k).

If scale is not 0, then the result is normalized such that trace(identity)=scale.

source
PauliStrings.trace_product_zFunction
trace_product_z(o1::AbstractOperator, o2::AbstractOperator; scale=0)

Efficiently compute <0|o1*o2|0>. If scale is not 0, then the result is normalized such that trace(identity) = scale.

source
PauliStrings.momentsFunction
moments(H::AbstractOperator, kmax::Int; start=1, scale=0)

Compute the first kmax moments of H. start is the first moment to start from.

If scale is not 0, then the result is normalized such that trace(identity)=scale.

source

Random operators

Construction

Base.:+Method
Base.:+(o::Operator, args::Tuple{Number, Vararg{Any}})
Base.:+(o::Operator, args::Tuple{Vararg{Any}})
Base.:+(o::Operator, term::Tuple{Number, String})
Base.:+(o::Operator, term::String)

Main functions to contruct spin operators. Identical signatures are available for -.

Examples

k-local terms can be added by adding a tuple to the operator. The first element of the tuple is an optional coefficient. The other element are couples (symbol,site) where symbol can be "X", "Y", "Z", "Sx", "Sy", "Sz", "S+", "S-" and site is an integer specifying the site on wich the symbol is acting.

A = Operator(4)
A += 2, "X",1,"X",2
A += 3, "Y",1,"X",2
A += "X",3,"X",4
A += 4,"Z",3
A += 5.2,"X",1,"Y",2,"Z",3
julia> A
(4.0 + 0.0im) 11Z1
(3.0 - 0.0im) YX11
(1.0 + 0.0im) 11XX
(2.0 + 0.0im) XX11
(5.2 - 0.0im) XYZ1

Full strings can also be added:

A = Operator(4)
A += 2, "1XXY"
A += 2im, "11Z1"
julia> A
(0.0 + 2.0im) 11Z1
(2.0 - 0.0im) 1XXY
source
PauliStrings.complete_basisFunction
complete_basis(N::Int)

Return a vector of all the strings supported on N spins.

source
complete_basis(N::Int, support::Vector{Int})

Return a vector of all the strings supported on N spins with support in support. The support is a vector of integers between 1 and N, indicating the sites where the string can be non-identity.

source
PauliStrings.k_local_basisFunction
k_local_basis(N::Int, k::Int; atmost=false)

Return a vector of all the k-local strings supported on N spins.

Example

julia> k_local_basis(3,2)
27-element Vector{PauliString}:
 XX1
 YX1
 ZX1
 XY1
 YY1
 ⋮
 1YY
 1ZY
 1XZ
 1YZ
 1ZZ
source
k_local_basis(N::Int, k::Int, support::Vector{Int}; atmost=false)

Return a vector of all the k-local strings supported on N spins with support in support. The support is a vector of integers between 1 and N, indicating the sites where the string can be non-identity.

source
PauliStrings.k_local_basis_1dFunction
k_local_basis_1d(N::Int, k::Int; translational_symmetry::Bool=false)

Generates the basis of all k-site Pauli strings on N qubits, build from X,Y,Z and identity.

Arguments

  • N::Int: Number of qubits
  • k::Int: Maximum support size (number of sites)
  • translational_symmetry::Bool=false: If true, return only unit cell operators as OperatorTS1D; otherwise return all translations as Operators

Returns

  • Vector{<:AbstractOperator}: Basis of k-site Pauli strings
source
PauliStrings.x_basisFunction
x_basis(N::Int)

Return a vector of all the strings supported on N spins with only x and identity.

Example

julia> x_basis(2)
4-element Vector{PauliString}:
 11
 X1
 1X
 XX
source
x_basis(N::Int, support::Vector{Int})

Return a vector of all the strings supported on N spins with only x and identity, and with support in support. The support is a vector of integers between 1 and N, indicating the sites where the string can be non-identity.

source
PauliStrings.y_basisFunction
y_basis(N::Int)

Return a vector of all the strings supported on N spins with only y and identity.

Example

julia> y_basis(2)
4-element Vector{PauliString}:
 11
 Y1
 1Y
 YY
source
y_basis(N::Int, support::Vector{Int})

Return a vector of all the strings supported on N spins with only y and identity, and with support in support. The support is a vector of integers between 1 and N, indicating the sites where the string can be non-identity.

source
PauliStrings.z_basisFunction
z_basis(N::Int)

Return a vector of all the strings supported on N spins with only z and identity.

Example

julia> z_basis(2)
4-element Vector{PauliString}:
 11
 Z1
 1Z
 ZZ
source
z_basis(N::Int, support::Vector{Int})

Return a vector of all the strings supported on N spins with only z and identity, and with support in support. The support is a vector of integers between 1 and N, indicating the sites where the string can be non-identity.

source
PauliStrings.majoranaFunction
majorana(N::Int, k::Int)

Return the k-th Majorana operator on N spins. There are 2N Majoranas supported on N spins. They all anticomute :

\[\{\gamma_i, \gamma_j\} = 2\delta_{ij}\]

Example

julia> majorana(4,1)
(1.0 + 0.0im) X111

julia> majorana(4,2)
(1.0 - 0.0im) Y111

julia> majorana(4,3)
(1.0 + 0.0im) ZX11

julia> majorana(4,4)
(1.0 - 0.0im) ZY11
source
PauliStrings.string_2dFunction
string_2d(args::Tuple{Vararg{Any}}, L1::Int, L2::Int; pbc=false)

Helper functions to construct 2d pauli strings. The tuple is composed of triplets the form (P,x,y,...) where P is one of "X", "Y", "Z", "Sx", "Sy", "Sz", "S-", "S+", and x, y label the position in the lattice. If pbc = true, the x and y coordinates will always be brough back to the range $[1, L_1]$ and $[1, L_2]$ respectively.

Example:

L1 = L2 = 2
A = Operator(L1 * L2)
A += 0.5 * string_2d(("Z", 1, 1, "Z", 2, 1), L1, L2) # Horizontal interaction
A += 1.0 * string_2d(("Z", 1, 1, "Z", 1, 2), L1, L2) # Vertical interaction
julia> A
(1.0 + 0.0im) Z1Z1
(0.5 + 0.0im) ZZ11
source

Tranlation symmetry

PauliStrings.PauliStringTSType
PauliStringTS{Ls, Ps, T<:Unsigned} <: AbstractPauliString

Type representing a translation symmetric sum of a Pauli string. The tuple Ls specifies the period in each dimension. The tuple Ps specifies which dimensions are periodic (true) or open (false). By default all dimensions are periodic.

source
PauliStrings.OperatorTSType
OperatorTS{Ls}(o::Operator)
OperatorTS{Ls,Ps}(o::Operator)

Construct an $n$-dimensional translation symmetric operator from o where Ls is a tuple of integers (L1, L2, ...) and Ps is an optional tuple of Bool values indicating which dimensions are periodic (default: all true). The resulting operator is equivalent to

\[O_\mathrm{TS} = \sum_T T^\dag O T\]

where $T$ are all translations on the L1×L2×… hypercube. So if you feed it an operator that is already a sum, you should afterwards normalize it by the number of sites.

To get a dense operator from this lazy sum representation, see resum. To get a single term, see representative.

source
PauliStrings.OperatorTS1DFunction

Initialize a zero 1D translation-invariant operator on N qubits.

source
OperatorTS1D(o::Operator; full=true)

Initialize a 1D translation invariant operator from an Operator $O=\sum_i o_i O_i$ where $O_i=T_i(O_0)$ and $T_i$ is the i-sites translation operator. Set full=true if passing $O$, an Operator that is supported on the whole chain (i.e converting from a translation symmetric Operator) Set full=false if passing $O_0$,a local term o such that the full operator is $O=\sum_i o_i T_i(O_0)$

source
PauliStrings.periodicflagsFunction
periodicflags(::Type{<:PauliStringTS})
periodicflags(p::PauliStringTS)

Get the tuple of periodic flags of a given PauliStringTS.

source
PauliStrings.representativeFunction
representative(p::PauliStringTS)

Returns a unique representative string of the translation symmetric sum of the Pauli string p.

source
representative(o::OperatorTS)

Returns a unique term of the symmetric sum represented by o.

source
PauliStrings.resumFunction
resum(o::OperatorTS)

Perform the symmetric sum represented by o to yield a dense Operator containing unsymmetrized PauliStrings.

source
PauliStrings.is_tsFunction
is_ts(o::Operator)

return true if o is translation symmetric in one dimension

source
is_ts(o::Operator, Ls::Tuple)
is_ts(o::Operator, Ls::Tuple, Ps::Tuple)

return true if o is translation symmetric on a hypercube with side lengths Ls. Ps is an optional tuple of Bool values indicating which dimensions are periodic (default: all true).

source
PauliStrings.is_ts2dFunction
is_ts2d(o::Operator, L1)
is_ts2d(o::Operator, L1, Ps::NTuple{2,Bool})

return true if o is translation symmetric on a rectangle with sidelengths L1 × qubitlength(o)÷L1. Ps is an optional tuple of Bool values indicating which dimensions are periodic (default: both true).

source

States

PauliStrings.expectMethod
expect(o::Operator, state::String)

Computes the expectation value <state|o|state>. State is a single binary string that represents a pure state in the computational basis.

source
PauliStrings.expectMethod
expect(o::Operator, in_state::String, out_state::String)

Computes the expectation value <out_state|o|in_state>. in_state and out_state are single binary strings that represents pure states in the computational basis.

source
PauliStrings.expect_productMethod
expect_product(o1::Operator, o2::Operator, state::String)

Computes the expectation value <state|o1*o2|state>. State is a single binary string that represents a pure state in the computational basis.

source

Circuits

PauliStrings.Circuits.CircuitType
Circuit(N::Int; max_strings=2^30, noise_amplitude=0)

Creates an empty quantum circuit with N qubits. max_strings is the maximum number of strings to keep in the operator. noise_amplitude is the amplitude of the noise gate.

source
Base.push!Method
push!(c::Circuit, gate::String, sites::Real...)

Adds a gate to the circuit c. The gate is specified by a string gate and a list of sites sites. The gates have the same naming convention as in Qiskit. Allowed gates are: "X", "Y", "Z", "H", "S", "T", "Tdg", "Phase", "CNOT", "Swap", "CX", "CY", "CZ", "CCX", "CSX", "CSXdg", "MCZ", "Noise".

source
Base.pushfirst!Method
pushfirst!(c::Circuit, gate::String, sites::Real...)

Adds a gate to the beginning of the circuit c.

source
PauliStrings.Circuits.XGateFunction
XGate(N::Int, i::Int)
YGate(N::Int, i::Int)
ZGate(N::Int, i::Int)
HGate(N::Int, i::Int)
SGate(N::Int, i::Int)
TGate(N::Int, i::Int)
TdgGate(N::Int, i::Int)
SXGate(N::Int, i::Int)

Creates a single qubit gate acting on qubit i of a N qubit system.

source
PauliStrings.Circuits.UGateFunction
UGate(N::Int, i::Int, theta::Real, phi::Real, lam::Real)

General 1-qubit rotation of qubit i of a N qubit system with Euler angles theta, phi, lam of form

$U(\theta, \phi, \lambda) =\begin{pmatrix}\cos\left(\theta/2\right) & -e^{i\lambda}\sin\left(\theta/2\right) \\e^{i\phi}\sin\left(\theta/2\right) & e^{i(\phi+\lambda)}\cos\left(\theta/2\right)\end{pmatrix}$

source
PauliStrings.Circuits.RXGateFunction
RXGate(N::Int, i::Int, phi::Real)
RYGate(N::Int, i::Int, theta::Real)
RZGate(N::Int, i::Int, phi::Real)

1-qubit rotation of qubit i of a N qubit system around specific axis.

source
PauliStrings.Circuits.CXGateFunction
CXGate(N::Int, i::Int, j::Int)
CYGate(N::Int, i::Int, j::Int)
CZGate(N::Int, i::Int, j::Int)
CNOTGate(N::Int, i::Int, j::Int)

Creates a two qubit gate with control qubit i and target qubit j of a N qubit system.

source
PauliStrings.Circuits.compileFunction
compile(c::Circuit)

Compiles the quantum circuit c into a unitary operator. Applies the gates in the order they were added. Applies noise gates if present and trim the operator to c.max_strings strings at each step.

source
PauliStrings.expectMethod
expect(c::Circuit, state::String)

Computes the expectation value <state|c|0>. State is a single binary string that represents a pure state in the computational basis.

source
PauliStrings.expectMethod
expect(c::Circuit, in_state::String, out_state::String)

Computes the expectation value of the state out_state after applying the circuit c to the state in_state.

source

I/O and conversion

PauliStrings.op_to_stringsMethod
op_to_strings(o::Operator)

takes an operator, return (coefs, strings) where coefs is a list of numbers and strings is a list of pauli string coefs[i] multiply strings[i]

source
Base.MatrixMethod
Matrix(o::Operator)

Convert an operator to a dense matrix.

source
PauliStrings.get_pauliFunction
get_pauli(o::Operator, i::Int)

Return an operator that represent the i-th pauli string of `o'. Does not return the string multiplied by the coefficient. Only the string.

source
Base.stringMethod
Base.string(x::PauliString)

Convert a PauliString to its string representation.

source

Symbolics

PauliStrings.simplify_operatorMethod
simplify_operator(o::Operator{P,Complex{Num}}) where {P}

Simplifies an Operator defined with symbolic coefficients. Uses Symbolics.simplify to simplify the symbolic expressions in each of the coefficients of o. Returns a new Operator.

source
PauliStrings.substitute_operatorMethod
substitute_operator(o::Operator{P,Complex{Num}}, dict::Dict) where {P}

Substitutes some or all of the variables in o according to the rule(s) in dict. If all the substitutions are to concrete numeric values, then it will return an Operator with Complex64 coefficients.

source

Other tools

PauliStrings.supportFunction
 support(p::PauliString)

Return the support of a Pauli string, i.e. the set of qubits on which it acts non-trivially.

source

Index