Documentation
Basics
PauliStrings.Operator — MethodOperator(N::Integer)Initialize a zero operator on N qubits.
PauliStrings.OperatorTS — MethodOperatorTS{Ls}(o::Operator)Construct an $n$-dimensional translation symmetric operator from o where Ls is a tuple of integers (L1, L2, ...) 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.
PauliStrings.qubitlength — Functionqubitlength(x::AbstractOperator)
qubitlength(::Type{<:AbstractOperator})Returns the number of qubits the operator acts on.
Truncation and noise
PauliStrings.add_noise — Methodadd_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
PauliStrings.add_noise — Methodadd_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.
PauliStrings.add_dephasing_noise — Methodadd_dephasing_noise(o::AbstractOperator, g::Real)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
PauliStrings.add_dephasing_noise — Methodadd_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$.
Base.truncate — Methodtruncate(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",4julia> A
(1.0 + 0.0im) ZZ1Z
(1.0 + 0.0im) XX11
julia> ps.truncate(A,2)
(1.0 + 0.0im) XX11PauliStrings.k_local_part — Methodk_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) XX11PauliStrings.trim — Methodtrim(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) XX11PauliStrings.prune — Method 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
PauliStrings.cutoff — Methodcutoff(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) ZZXXAlgorithms
PauliStrings.rk4 — Methodrk4(H::AbstractOperator, O::AbstractOperator, dt::Real; hbar::Real=1, heisenberg=true, M=2^20, keep::Operator=Operator(0))Single step of Runge–Kutta-4 with time independant Hamiltonian. Returns O(t+dt). Set heisenberg=true for evolving an observable in the heisenberg picture. If heisenberg=false then it is assumed that O is a density matrix. M is the number of strings to keep.
PauliStrings.rk4 — Methodrk4(H::Function, O::AbstractOperator, dt::Real, t::Real; hbar::Real=1, heisenberg=true, M=2^20, keep::Operator=Operator(0))Single step of Runge–Kutta-4 with time dependant Hamiltonian. Returns O(t+dt). H is a function that takes a number (time) and returns an operator. Set heisenberg=true for evolving an observable in the heisenberg picture. If heisenberg=false then it is assumed that O is a density matrix. M is the number of strings to keep.
PauliStrings.rk4_lindblad — Methodrk4_lindblad(H::AbstractOperator, O::AbstractOperator, dt::Real, L; hbar::Real=1, heisenberg=true, M=2^20, keep::Operator=Operator(0), 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. Set heisenberg=true for evolving an observable in the heisenberg picture. If heisenberg=false then it is assumed that O is a density matrix.
PauliStrings.lanczos — Methodlanczos(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
Operations
Base.:+ — MethodBase.:+(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", 3julia> 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) 1111Base.:- — MethodBase.:-(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
Base.:* — MethodBase.:*(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", 3julia> 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) XYZ1Base.:^ — MethodBase.:^(o::Operator, k::Int)kth power of o.
PauliStrings.commutator — Methodcommutator(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) Y111PauliStrings.anticommutator — Methodanticommutator(o1::Operator, o2::Operator; kwargs...)Commutator of two operators. This is faster than doing o1*o2 + o2*o1.
Base.:/ — MethodBase.:/(o::AbstractOperator, a::Number)Divide an operator by a number
PauliStrings.compress — Methodcompress(o::AbstractOperator)Accumulate repeated terms
PauliStrings.trace — Methodtrace(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.0imPauliStrings.diag — Methoddiag(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) Z11ZPauliStrings.opnorm — Methodopnorm(o::AbstractOperator; normalize=false)Frobenius norm. If normalize is true, return the trace divided by sqrt(2^N).
Example
julia> A = Operator(4)
julia> A += 2,"X",2
julia> A += 1,"Z",1,"Z",3
julia> opnorm(A)
8.94427190999916PauliStrings.dagger — Methoddagger(o::AbstractOperator)Conjugate transpose. ' also works.
Example
A = Operator(3)
A += 1im,"X",2
A += 1,"Z",1,"Z",3julia> A
(1.0 + 0.0im) Z1Z
(0.0 + 1.0im) 1X1
julia> dagger(A)
(1.0 - 0.0im) Z1Z
(0.0 - 1.0im) 1X1
julia> A'
(1.0 - 0.0im) Z1Z
(0.0 - 1.0im) 1X1PauliStrings.ptrace — Methodptrace(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) XY11ZPower and moments
Base.:^ — MethodBase.:^(o::Operator, k::Int)kth power of o.
PauliStrings.trace_product — Methodtrace_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.
PauliStrings.trace_product — Methodtrace_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.
PauliStrings.trace_product — Methodtrace_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.
PauliStrings.trace_product_z — Methodtrace_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.
PauliStrings.moments — Methodmoments(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.
Random operators
PauliStrings.rand_local1 — Methodrand_local2(N::Int)Random 1-local operator
PauliStrings.rand_local2 — Methodrand_local2(N::Int)Random 2-local operator
PauliStrings.rand_local1_TS1D — Methodrand_local1_TS1D(N::Int)Random 1-local OperatorTS in one dimension
PauliStrings.rand_local2_TS1D — Methodrand_local2_TS1D(N::Int)Random 2-local OperatorTS in one dimension
Construction
Base.:+ — MethodBase.:+(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",3julia> 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) XYZ1Full 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) 1XXYPauliStrings.all_strings — Methodall_strings(N::Int)Return the sum of all the strings supported on N spins, with coeficients 1
PauliStrings.all_k_local — Methodall_k_local(N::Int, k::Int; atmost=false)Return the sum of all the k-local strings supported on N spins, with coeficients 1. These are k-local only strings, not including strings shorter than k.
Example
julia> all_k_local(2, 1)
(1.0 + 0.0im) X1
(1.0 + 0.0im) 1X
(1.0 + 0.0im) Z1
(1.0 - 0.0im) Y1
(1.0 + 0.0im) 1Z
(1.0 - 0.0im) 1YPauliStrings.all_x — Methodall_x(N::Int)Return the sum of all the strings supported on N spins with only x and with coeficients 1
Example
julia> all_x(2)
(1.0 + 0.0im) 11
(1.0 + 0.0im) X1
(1.0 + 0.0im) 1X
(1.0 + 0.0im) XXPauliStrings.all_y — Methodall_y(N::Int)Return the sum of all the strings supported on N spins with only y and with coeficients 1
Example
julia> all_y(2)
(1.0 + 0.0im) 11
(1.0 - 0.0im) Y1
(1.0 - 0.0im) 1Y
(1.0 - 0.0im) YYPauliStrings.all_z — Methodall_z(N::Int)Return the sum of all the strings supported on N spins with only z and with coeficients 1
Example
julia> all_z(2)
(1.0 + 0.0im) 11
(1.0 + 0.0im) Z1
(1.0 + 0.0im) 1Z
(1.0 + 0.0im) ZZPauliStrings.majorana — Methodmajorana(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) ZY11PauliStrings.string_2d — Methodstring_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 interactionjulia> A
(1.0 + 0.0im) Z1Z1
(0.5 + 0.0im) ZZ11Tranlation symmetry
PauliStrings.PauliStringTS — MethodPauliStringTS{Ls}(p::PauliString)Construct a translation symmetric sum of a Pauli string p. Ls is a tuple that specifies the periodicity in each dimension.
PauliStrings.OperatorTS — MethodOperatorTS{Ls}(o::Operator)Construct an $n$-dimensional translation symmetric operator from o where Ls is a tuple of integers (L1, L2, ...) 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.
PauliStrings.OperatorTS1D — MethodOperatorTS1D(N::Integer)Initialize a zero 1D translation-invariant operator on N qubits.
PauliStrings.OperatorTS1D — MethodOperatorTS1D(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)$
PauliStrings.representative — Methodrepresentative(o::OperatorTS)Returns a unique term of the symmetric sum represented by o.
PauliStrings.representative — Methodrepresentative(p::PauliStringTS)Returns a unique representative string of the translation symmetric sum of the Pauli string p.
PauliStrings.resum — Methodresum(o::OperatorTS)Perform the symmetric sum represented by o to yield a dense Operator containing unsymmetrized PauliStrings.
PauliStrings.is_ts — Methodis_ts(o::Operator)return true if o is translation symmetric in one dimension
PauliStrings.is_ts — Methodis_ts(o::Operator, Ls::Tuple)return true if o is translation symmetric on a hypercube with side lengths Ls.
PauliStrings.is_ts2d — Methodis_ts2d(o::Operator, L1)return true if o is translation symmetric on a rectangle with sidelengths L1 × qubitlength(o)÷L1.
States
PauliStrings.trace_zpart — Methodtrace_zpart(o::Operator)Computes <0|o|0>.
PauliStrings.expect — Methodexpect(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.
PauliStrings.expect — Methodexpect(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.
PauliStrings.expect_product — Methodexpect_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.
Circuits
PauliStrings.Circuits.Circuit — MethodCircuit(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.
Base.push! — Methodpush!(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".
Base.pushfirst! — Methodpushfirst!(c::Circuit, gate::String, sites::Real...)Adds a gate to the beginning of the circuit c.
PauliStrings.Circuits.XGate — MethodXGate(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.
PauliStrings.Circuits.UGate — MethodUGate(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}$
PauliStrings.Circuits.RXGate — MethodRXGate(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.
PauliStrings.Circuits.PhaseGate — MethodPhaseGate(N::Int, i::Int, theta::Real)Creates a phase gate acting on qubit i of a N qubit system with phase theta.
PauliStrings.Circuits.CPhaseGate — MethodCPhaseGate(N::Int, i::Int, j::Int, theta::Real)Controled phase gate with control qubit i and target qubit j of a N qubit system.
PauliStrings.Circuits.CXGate — MethodCXGate(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.
PauliStrings.Circuits.SwapGate — MethodSwapGate(N::Int, i::Int, j::Int)Creates a swap gate between qubits i and j of a N qubit system.
PauliStrings.Circuits.CSXGate — MethodCSXGate(N::Int, i::Int, j::Int)
CSXdgGate(N::Int, i::Int, j::Int)Controlled sqrt X gate and its dagger
PauliStrings.Circuits.CCXGate — MethodCCXGate(N::Int, i::Int, j::Int, k::Int)Tofolli gate with control qubits i and j and target qubit k of a N qubit system.
PauliStrings.Circuits.MCZGate — MethodMCZGate(N::Int, sites::Int...)Creates a multi-controlled Z gate acting on sites qubits of a N qubit system.
PauliStrings.Circuits.XXPlusYYGate — MethodXXPlusYYGate(N::Int, i::Int, j::Int, theta::Real, beta::Real)XX+YY gate between qubits i and j of a N qubit system.
PauliStrings.Circuits.grover_diffusion — Methodgrover_diffusion(N::Int, sites::Int...)Creates the Grover diffusion operator acting on sites qubits of a N qubit system.
PauliStrings.Circuits.compile — Methodcompile(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.
PauliStrings.expect — Methodexpect(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.
PauliStrings.expect — Methodexpect(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.
I/O and conversion
PauliStrings.op_to_strings — Methodop_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]
PauliStrings.get_coeffs — Methodget_coeffs(o::AbstractOperator)Return the list of coefficient in front of each strings.
PauliStrings.set_coeffs — Methodset_coeffs(o::Operator, coefs::Vector{T}) where T <: NumberSets the coefficient of o to coefs. Inplace.
A = Operator(4)
A += 2, "1XXY"
A += 3, "11Z1"julia> A
(3.0 + 0.0im) 11Z1
(2.0 - 0.0im) 1XXY
julia> set_coeffs(A, [5,6])
julia> A
(5.0 + 0.0im) 11Z1
(6.0 - 0.0im) 1XXYPauliStrings.op_to_dense — Methodop_to_dense(o::Operator)Convert an operator to a dense matrix.
Base.Matrix — MethodMatrix(o::Operator)Convert an operator to a dense matrix.
PauliStrings.get_coeff — Methodget_coeff(o::Operator{P}, p::P) where {P}Return the coefficient of the string p in o.
PauliStrings.get_pauli — Methodget_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.
Other tools
PauliStrings.compress — Methodcompress(o::AbstractOperator)Accumulate repeated terms
PauliStrings.xcount — Methodxcount(p::PauliString)Count the number of X operators in a string.
PauliStrings.ycount — Methodycount(p::PauliString)Count the number of Y operators in a string.
PauliStrings.zcount — Methodzcount(p::PauliString)Count the number of Z operators in a string.
Index
Base.MatrixPauliStrings.Circuits.CircuitPauliStrings.OperatorPauliStrings.OperatorTSPauliStrings.OperatorTSPauliStrings.PauliStringTSBase.:*Base.:+Base.:+Base.:-Base.:/Base.:^Base.:^Base.push!Base.pushfirst!Base.truncatePauliStrings.Circuits.CCXGatePauliStrings.Circuits.CPhaseGatePauliStrings.Circuits.CSXGatePauliStrings.Circuits.CXGatePauliStrings.Circuits.MCZGatePauliStrings.Circuits.PhaseGatePauliStrings.Circuits.RXGatePauliStrings.Circuits.SwapGatePauliStrings.Circuits.UGatePauliStrings.Circuits.XGatePauliStrings.Circuits.XXPlusYYGatePauliStrings.Circuits.compilePauliStrings.Circuits.grover_diffusionPauliStrings.OperatorTS1DPauliStrings.OperatorTS1DPauliStrings.add_dephasing_noisePauliStrings.add_dephasing_noisePauliStrings.add_noisePauliStrings.add_noisePauliStrings.all_k_localPauliStrings.all_stringsPauliStrings.all_xPauliStrings.all_yPauliStrings.all_zPauliStrings.anticommutatorPauliStrings.commutatorPauliStrings.compressPauliStrings.compressPauliStrings.cutoffPauliStrings.daggerPauliStrings.diagPauliStrings.expectPauliStrings.expectPauliStrings.expectPauliStrings.expectPauliStrings.expect_productPauliStrings.get_coeffPauliStrings.get_coeffsPauliStrings.get_pauliPauliStrings.is_tsPauliStrings.is_tsPauliStrings.is_ts2dPauliStrings.k_local_partPauliStrings.lanczosPauliStrings.majoranaPauliStrings.momentsPauliStrings.op_to_densePauliStrings.op_to_stringsPauliStrings.opnormPauliStrings.prunePauliStrings.ptracePauliStrings.qubitlengthPauliStrings.rand_local1PauliStrings.rand_local1_TS1DPauliStrings.rand_local2PauliStrings.rand_local2_TS1DPauliStrings.representativePauliStrings.representativePauliStrings.resumPauliStrings.rk4PauliStrings.rk4PauliStrings.rk4_lindbladPauliStrings.set_coeffsPauliStrings.string_2dPauliStrings.tracePauliStrings.trace_productPauliStrings.trace_productPauliStrings.trace_productPauliStrings.trace_product_zPauliStrings.trace_zpartPauliStrings.trimPauliStrings.xcountPauliStrings.ycountPauliStrings.zcount