Full API

States

Stabilizer states can be represented with the Stabilizer, Destabilizer, MixedStabilizer, and MixedDestabilizer tableau data structures. You probably want to use MixedDestabilizer which supports the widest set of operations.

Moreover, a MixedDestabilizer can be stored inside a Register together with a set of classical bits in which measurement results can be written.

Lastly, for Pauli frame simulations there is the PauliFrame type, a tableau in which each row represents a different Pauli frame.

There are convenience constructors for common types of states and operators.

Operations

Acting on quantum states can be performed either:

  • In a "linear algebra" language where unitaries, measurements, and other operations have separate interfaces. This is an explicitly deterministic lower-level interface, which provides a great deal of control over how tableaux are manipulated. See the Stabilizer Tableau Algebra Manual as a primer on these approaches.
  • Or in a "circuit" language, where the operators (and measurements and noise) are represented as circuit gates. This is a higher-level interface in which the outcome of an operation can be stochastic. The API for it is centered around the apply! function. Particularly useful for Monte Carlo simulations and Perturbative Expansion Symbolic Results.

See the full list of operations for a list of implemented operations.

Autogenerated API list

QuantumClifford.BellMeasurementType

A Bell measurement performing the correlation measurement corresponding to the given pauli projections on the qubits at the selected indices.

source
QuantumClifford.CliffordOperatorType

Clifford Operator specified by the mapping of the basis generators.

julia> tCNOT
X₁ ⟼ + XX
X₂ ⟼ + _X
Z₁ ⟼ + Z_
Z₂ ⟼ + ZZ

julia> phase_gate = C"Y
                      Z"
X₁ ⟼ + Y
Z₁ ⟼ + Z

julia> stab = S"XI
                IZ";


julia> entangled = tCNOT*stab
+ XX
+ ZZ

julia> CliffordOperator(T"YY")
ERROR: DimensionMismatch: Input tableau should be of size 2n×n (top half is the X mappings and the bottom half are the Z mappings).
[...]

Destabilizer can also be converted.

julia> d = Destabilizer(S"Y")
𝒟ℯ𝓈𝓉𝒶𝒷
+ Z
𝒮𝓉𝒶𝒷
+ Y

julia> CliffordOperator(d)
X₁ ⟼ + Z
Z₁ ⟼ + Y
source
QuantumClifford.DestabilizerType

A tableau representation of a pure stabilizer state. The tableau tracks the destabilizers as well, for efficient projections. On initialization there are no checks that the provided state is indeed pure. This enables the use of this data structure for mixed stabilizer state, but a better choice would be to use MixedDestabilizer.

source
QuantumClifford.MixedDestabilizerType

A tableau representation for mixed stabilizer states that keeps track of the destabilizers in order to provide efficient projection operations.

The rank r of the n-qubit tableau is tracked, either so that it can be used to represent a mixed stabilizer state, or so that it can be used to represent an n-r logical-qubit code over n physical qubits. The "logical" operators are tracked as well.

When the constructor is called on an incomplete Stabilizer it automatically calculates the destabilizers and logical operators, following chapter 4 of (Gottesman, 1997). Under the hood the conversion uses the canonicalize_gott! canonicalization. That canonicalization permutes the columns of the tableau, but we automatically undo the column permutation in the preparation of a MixedDestabilizer so that qubits are not reindexed. The boolean keyword arguments undoperm and reportperm can be used to control this behavior and to report the permutations explicitly.

See also: stabilizerview, destabilizerview, logicalxview, logicalzview

source
QuantumClifford.PauliFrameType
struct PauliFrame{T, S} <: QuantumClifford.AbstractQCState

This is a wrapper around a tableau. This "frame" tableau is not to be viewed as a normal stabilizer tableau, although it does conjugate the same under Clifford operations. Each row in the tableau refers to a single frame. The row represents the Pauli operation by which the frame and the reference differ.

source
QuantumClifford.PauliFrameMethod
PauliFrame(
    frames,
    qubits,
    measurements
) -> PauliFrame{Stabilizer{QuantumClifford.Tableau{Vector{UInt8}, LinearAlgebra.Adjoint{UInt64, Matrix{UInt64}}}}}

Prepare an empty set of Pauli frames with the given number of frames and qubits. Preallocates spaces for measurement number of measurements.

source
QuantumClifford.PauliOperatorType

A multi-qubit Pauli operator ($±\{1,i\}\{I,Z,X,Y\}^{\otimes n}$).

A Pauli can be constructed with the P custom string macro or by building up one through products and tensor products of smaller operators.

julia> pauli3 = P"-iXYZ"
-iXYZ

julia> pauli4 = 1im * pauli3 ⊗ X
+ XYZX

julia> Z*X
+iY

We use a typical F(2,2) encoding internally. The X and Z bits are stored in a single concatenated padded array of UInt chunks of a bit array.

julia> p = P"-IZXY";


julia> p.xz
2-element Vector{UInt64}:
 0x000000000000000c
 0x000000000000000a

You can access the X and Z bits through getters and setters or through the xview, zview, xbit, and zbit functions.

julia> p = P"XYZ"; p[1]
(true, false)

julia> p[1] = (true, true); p
+ YYZ
source
QuantumClifford.RegisterType

A register, representing the state of a computer including both a tableaux and an array of classical bits (e.g. for storing measurement results)

source
QuantumClifford.ResetType

Reset the specified qubits to the given state.

Be careful, this operation implies first tracing out the qubits, which can lead to mixed states if these qubits were entangled with the rest of the system.

See also: sMRZ

source
QuantumClifford.SingleQubitOperatorType

A "symbolic" general single-qubit operator which permits faster multiplication than an operator expressed as an explicit tableau.

julia> op = SingleQubitOperator(2, true, true, true, false, true, true) # Tableau components and phases
SingleQubitOperator on qubit 2
X₁ ⟼ - Y
Z₁ ⟼ - X

julia> typeof(op)
SingleQubitOperator

julia> t_op = CliffordOperator(op, 3) # Transforming it back into an explicit tableau representation (specifying the size)
X₁ ⟼ + X__
X₂ ⟼ - _Y_
X₃ ⟼ + __X
Z₁ ⟼ + Z__
Z₂ ⟼ - _X_
Z₃ ⟼ + __Z

julia> typeof(t_op)
CliffordOperator{QuantumClifford.Tableau{Vector{UInt8}, Matrix{UInt64}}}

julia> CliffordOperator(op, 1, compact=true) # You can also extract just the non-trivial part of the tableau
X₁ ⟼ - Y
Z₁ ⟼ - X

See also: sHadamard, sPhase, sId1, sX, sY, sZ, CliffordOperator

Or simply consult subtypes(QuantumClifford.AbstractSingleQubitOperator) and subtypes(QuantumClifford.AbstractTwoQubitOperator) for a full list. You can think of the s prefix as "symbolic" or "sparse".

source
QuantumClifford.SparseGateType

A Clifford gate, applying the given cliff operator to the qubits at the selected indices.

apply!(state, cliff, indices) and apply!(state, SparseGate(cliff, indices)) give the same result.

source
QuantumClifford.StabMixtureType
mutable struct StabMixture{T, F}

Represents mixture ∑ ϕᵢⱼ Pᵢ ρ Pⱼ† where ρ is a pure stabilizer state.

julia> StabMixture(S"-X")
A mixture ∑ ϕᵢⱼ Pᵢ ρ Pⱼ† where ρ is
𝒟ℯ𝓈𝓉𝒶𝒷
+ Z
𝒮𝓉𝒶𝒷
- X
with ϕᵢⱼ | Pᵢ | Pⱼ:
 1.0+0.0im | + _ | + _

julia> pcT
A unitary Pauli channel P = ∑ ϕᵢ Pᵢ with the following branches:
with ϕᵢ | Pᵢ
 0.853553+0.353553im | + _
 0.146447-0.353553im | + Z

julia> apply!(StabMixture(S"-X"), pcT)
A mixture ∑ ϕᵢⱼ Pᵢ ρ Pⱼ† where ρ is
𝒟ℯ𝓈𝓉𝒶𝒷
+ Z
𝒮𝓉𝒶𝒷
- X
with ϕᵢⱼ | Pᵢ | Pⱼ:
 0.0+0.353553im | + _ | + Z
 0.0-0.353553im | + Z | + _
 0.853553+0.0im | + _ | + _
 0.146447+0.0im | + Z | + Z

See also: PauliChannel

source
QuantumClifford.StabilizerType

Stabilizer, i.e. a list of commuting multi-qubit Hermitian Pauli operators.

Instances can be created with the S custom string macro or as direct sum of other stabilizers.

Stabilizers and Destabilizers

In many cases you probably would prefer to use the MixedDestabilizer data structure, as it caries a lot of useful additional information, like tracking rank and destabilizer operators. Stabilizer has mostly a pedagogical value, and it is also used for slightly faster simulation of a particular subset of Clifford operations.

julia> s = S"XXX
             ZZI
             IZZ"
+ XXX
+ ZZ_
+ _ZZ

julia> s⊗s
+ XXX___
+ ZZ____
+ _ZZ___
+ ___XXX
+ ___ZZ_
+ ____ZZ

It has an indexing API, looking like a list of PauliOperators.

julia> s[2]
+ ZZ_

Pauli operators can act directly on the a stabilizer.

julia> P"YYY" * s
- XXX
+ ZZ_
+ _ZZ

There are a number of ways to create a Stabilizer, including:

  • generate Stabilizers from a list of Pauli operators
julia> Stabilizer([P"XX", P"ZZ"])
+ XX
+ ZZ
  • generate Stabilizers from boolean matrices
julia> a = [true true; false false]; b = [false true; true true];

julia> Stabilizer(a, b)
+ XY
+ ZZ

julia> Stabilizer([0x0, 0x2], a, b)
+ XY
- ZZ
  • initialize an empty Stabilizer and fill it through indexing
julia> s = zero(Stabilizer, 2)
+ __
+ __

julia> s[1,1] = (true, false); s
+ X_
+ __

There are no automatic checks for correctness (i.e. independence of all rows, commutativity of all rows, hermiticity of all rows). The rank (number of rows) is permitted to be less than the number of qubits (number of columns): canonilization, projection, etc. continue working in that case. To great extent this library uses the Stabilizer data structure simply as a tableau. This might be properly abstracted away in future versions.

See also: PauliOperator, canonicalize!

source
QuantumClifford.UnitaryPauliChannelType

A Pauli channel datastructure, mainly for use with StabMixture.

More convenient to use than PauliChannel when you know your Pauli channel is unitary.

julia> Tgate = UnitaryPauliChannel(
           (I, Z),
           ((1+exp(im*π/4))/2, (1-exp(im*π/4))/2)
       )
A unitary Pauli channel P = ∑ ϕᵢ Pᵢ with the following branches:
with ϕᵢ | Pᵢ
 0.853553+0.353553im | + _
 0.146447-0.353553im | + Z

julia> PauliChannel(Tgate)
Pauli channel ρ ↦ ∑ ϕᵢⱼ Pᵢ ρ Pⱼ† with the following branches:
with ϕᵢⱼ | Pᵢ | Pⱼ:
 0.853553+0.0im | + _ | + _
 0.0+0.353553im | + _ | + Z
 0.0-0.353553im | + Z | + _
 0.146447+0.0im | + Z | + Z
source
QuantumClifford.VerifyOpType

A "probe" to verify that the state of the qubits corresponds to a desired good_state, e.g. at the end of the execution of a circuit.

source
QuantumClifford.sMRZType

Measure a qubit in the Z basis and reset to the |0⟩ state.

It does not trace out the qubit!

As described below there is a difference between measuring the qubit (followed by setting it to a given known state) and "tracing out" the qubit. By reset here we mean "measuring and setting to a known state", not "tracing out".

julia> s = MixedDestabilizer(S"XXX ZZI IZZ") # |000⟩+|111⟩
𝒟ℯ𝓈𝓉𝒶𝒷
+ Z__
+ _X_
+ __X
𝒮𝓉𝒶𝒷━
+ XXX
+ ZZ_
+ Z_Z

julia> traceout!(copy(s), 1) # = I⊗(|00⟩⟨00| + |11⟩⟨11|)
𝒟ℯ𝓈𝓉𝒶𝒷
+ _X_
𝒳ₗ━━━
+ _XX
+ Z__
𝒮𝓉𝒶𝒷━
+ _ZZ
𝒵ₗ━━━
+ Z_Z
+ XXX

julia> projectZ!(traceout!(copy(s), 1), 1)[1] # = |000⟩⟨000|+|011⟩⟨011| or |100⟩⟨100|+|111⟩⟨111| (use projectZrand! to actually get a random result)
𝒟ℯ𝓈𝓉𝒶𝒷
+ _X_
+ XXX
𝒳ₗ━━━
+ _XX
𝒮𝓉𝒶𝒷━
+ _ZZ
+ Z__
𝒵ₗ━━━
+ Z_Z

julia> projectZ!(copy(s), 1)[1] # = |000⟩ or |111⟩ (use projectZrand! to actually get a random result)
𝒟ℯ𝓈𝓉𝒶𝒷
+ XXX
+ _X_
+ __X
𝒮𝓉𝒶𝒷━
+ Z__
+ ZZ_
+ Z_Z
julia> apply!(Register(copy(s)), sMRZ(1)) |> quantumstate # |000⟩ or |011⟩, depending on randomization
𝒟ℯ𝓈𝓉𝒶𝒷
+ XXX
+ _X_
+ __X
𝒮𝓉𝒶𝒷━
+ Z__
- ZZ_
- Z_Z

See also: Reset, sMZ

source
QuantumClifford.PauliErrorFunction

A convenient constructor for various types of Pauli errors, that can be used as circuit gates in simulations. Returns more specific types when necessary.

source
QuantumClifford.PauliErrorMethod

"Construct a gate operation that applies a biased Pauli error on all qubits independently, each with probabilities px, py, pz. Note that the probability of any error occurring is px+py+pz. Because of this, PauliError(1, p) is equivalent to PauliError(1,p/3,p/3,p/3). Similarly, if one wanted to exclude Z errors from PauliError(1,p/3,p/3,p/3) while mainting the same rate of X errors, one could write PauliError(1, p*2/3, 0, 0) (in the sense that Y errors can be interpreted as an X and a Z happening at the same time).

source
QuantumClifford.PauliErrorMethod

"Construct a gate operation that applies a biased Pauli error on qubit q with independent probabilities px, py, pz. Note that the probability of any error occurring is px+py+pz. Because of this, PauliError(1, p) is equivalent to PauliError(1,p/3,p/3,p/3). Similarly, if one wanted to exclude Z errors from PauliError(1,p/3,p/3,p/3) while mainting the same rate of X errors, one could write PauliError(1, p*2/3, 0, 0) (in the sense that Y errors can be interpreted as an X and a Z happening at the same time).

source
QuantumClifford.applybranchesFunction

Compute all possible new states after the application of the given operator. Reports the probability of each one of them. Deterministic (as it reports all branches of potentially random processes), part of the Perturbative Expansion interface.

source
QuantumClifford.applynoise!Function

A method modifying a given state by applying the corresponding noise model. It is non-deterministic, part of the Noise interface.

source
QuantumClifford.bellFunction

Prepare one or more Bell pairs (with optional phases).

julia> bell()
+ XX
+ ZZ

julia> bell(2)
+ XX__
+ ZZ__
+ __XX
+ __ZZ

julia> bell((true, false))
- XX
+ ZZ

julia> bell([true, false, true, true])
- XX__
+ ZZ__
- __XX
- __ZZ
source
QuantumClifford.bigramMethod
bigram(
    state::QuantumClifford.AbstractStabilizer;
    clip
) -> Matrix{Int64}

Get the bigram of a tableau.

It is the list of endpoints of a tableau in the clipped gauge.

If clip=true (the default) the tableau is converted to the clipped gauge in-place before calculating the bigram. Otherwise, the clip gauge conversion is skipped (for cases where the input is already known to be in the correct gauge).

Introduced in (Nahum et al., 2017), with a more detailed explanation of the algorithm in (Li et al., 2019) and (Gullans et al., 2021).

See also: canonicalize_clip!

source
QuantumClifford.canonicalize!Method
canonicalize!(
    state::QuantumClifford.AbstractStabilizer;
    phases,
    ranks
) -> Union{Tuple{QuantumClifford.AbstractStabilizer, Int64, Int64}, QuantumClifford.AbstractStabilizer}

Canonicalize a stabilizer (in place).

Assumes the input is a valid stabilizer (all operators commute and have real phases). It permits redundant generators and identity generators.

julia> ghz = S"XXXX
               ZZII
               IZZI
               IIZZ";


julia> canonicalize!(ghz)
+ XXXX
+ Z__Z
+ _Z_Z
+ __ZZ

julia> canonicalize!(S"XXXX
                       IZZI
                       IIZZ")
+ XXXX
+ _Z_Z
+ __ZZ

Not all rows in the tableau in the next example are independent:

julia> canonicalize!(S"XXXX
                       ZZII
                       IZZI
                       IZIZ
                       IIZZ")
+ XXXX
+ Z__Z
+ _Z_Z
+ __ZZ
+ ____

In cases of lower rank, more advanced tableau structures might be better. For instance the MixedStabilizer or MixedDestabilizer structures (you can read more about them in the Data Structures section of the documentation).

If phases=false is set, the canonicalization does not track the phases in the tableau, leading to significant (constant factor) speedup.

julia> s = S"-ZX
              XZ"
- ZX
+ XZ

julia> canonicalize!(copy(s), phases=false)
- XZ
+ ZX

julia> canonicalize!(copy(s))
+ XZ
- ZX

If ranks=true is set, the last pivot indices for the X and Z stage of the canonicalization are returned as well.

julia> s = S"XXXX
             ZZII
             IZIZ
             ZIIZ";


julia> _, ix, iz = canonicalize!(s, ranks=true); ix, iz
(1, 3)

julia> s
+ XXXX
+ Z__Z
+ _Z_Z
+ ____

Based on (Garcia et al., 2012).

See also: canonicalize_rref!, canonicalize_gott!

source
QuantumClifford.canonicalize_clip!Method
canonicalize_clip!(
    state::QuantumClifford.AbstractStabilizer;
    phases
) -> QuantumClifford.AbstractStabilizer

Fix the clipped gauge of a stabilizer (in place).

Assumes the input is a valid full-rank stabilizer (all operators commute and have real phases).

julia> s = S"- X_ZX_X
             + XXYZ__
             - YZ_Z_X
             - XZX__Y
             + _Z_Y_Y
             - ____Z_";


julia> canonicalize_clip!(s)
- X_XY__
+ YZY___
+ _XZX__
- _ZYX_Z
- __YZ_X
- ____Z_

If phases=false is set, the canonicalization does not track the phases in the tableau, leading to a significant speedup.

Introduced in (Nahum et al., 2017), with a more detailed explanation of the algorithm in Appendix A of (Li et al., 2019)

See also: canonicalize!, canonicalize_rref!, canonicalize_gott!.

source
QuantumClifford.canonicalize_gott!Method

Inplace Gottesman canonicalization of a tableau.

This uses different canonical form from canonicalize!. It is used in the computation of the logical X and Z operators of a MixedDestabilizer.

It returns the (in place) modified state, the indices of the last pivot of both Gaussian elimination steps, and the permutations that have been used to put the X and Z tableaux in standard form.

Based on (Gottesman, 1997).

See also: canonicalize!, canonicalize_rref!

source
QuantumClifford.canonicalize_rref!Method
canonicalize_rref!(
    state::QuantumClifford.AbstractStabilizer,
    colindices;
    phases
) -> Tuple{QuantumClifford.AbstractStabilizer, Any}

Canonicalize a stabilizer (in place) along only some columns.

This uses different canonical form from canonicalize!. It also indexes in reverse in order to make its use in traceout! more efficient. Its use in traceout! is its main application.

It returns the (in place) modified state and the index of the last pivot.

Based on (Audenaert and Plenio, 2005).

See also: canonicalize!, canonicalize_gott!

source
QuantumClifford.centralizerMethod

For a given set of Paulis (in the form of a Tableau), return the subset of Paulis that commute with all Paulis in set.

julia> centralizer(T"XX ZZ _Z")
+ ZZ
source
QuantumClifford.clifford_cardinalityMethod

The size of the Clifford group 𝒞 over a given number of qubits, possibly modulo the phases.

For n qubits, not accounting for phases is 2ⁿⁿΠⱼ₌₁ⁿ(4ʲ-1). There are 4ⁿ different phase configurations.

julia> clifford_cardinality(7)
457620995529680351512370381586432000

When not accounting for phases (phases = false) the result is the same as the size of the Symplectic group Sp(2n) ≡ 𝒞ₙ/𝒫ₙ, where 𝒫ₙ is the Pauli group over n qubits.

julia> clifford_cardinality(7, phases=false)
27930968965434591767112450048000

See also: enumerate_cliffords.

source
QuantumClifford.commFunction

Check whether two operators commute.

0x0 if they commute, 0x1 if they anticommute.

julia> P"XX"*P"ZZ", P"ZZ"*P"XX"
(- YY, - YY)

julia> comm(P"ZZ", P"XX")
0x00

julia> comm(P"IZ", P"XX")
0x01

See also: comm!

source
QuantumClifford.compactify_circuitMethod

Convert a list of gates to a more optimized "sum type" format which permits faster dispatch.

Generally, this should be called on a circuit before it is used in a simulation.

source
QuantumClifford.contractorMethod

Return the subset of Paulis in a Stabilizer that have identity operators on all qubits corresponding to the given subset, without the entries corresponding to subset.

julia> contractor(S"_X X_", [1])
+ X
source
QuantumClifford.delete_columnsMethod

Return the given stabilizer without all the qubits in the given iterable.

The resulting tableaux is not guaranteed to be valid (to retain its commutation relationships).

julia> delete_columns(S"XYZ YZX ZXY", [1,3])
+ Y
+ Z
+ X

See also: traceout!

source
QuantumClifford.enumerate_single_qubit_gatesMethod

Generate a symbolic single-qubit gate given its index. Optionally, set non-trivial phases.

julia> enumerate_single_qubit_gates(6)
sPhase on qubit 1
X₁ ⟼ + Y
Z₁ ⟼ + Z

julia> enumerate_single_qubit_gates(6, qubit=2, phases=(true, true))
SingleQubitOperator on qubit 2
X₁ ⟼ - Y
Z₁ ⟼ - Z

See also: enumerate_cliffords.

source
QuantumClifford.fastcolumnFunction

Convert a tableau to a memory layout that is fast for column operations.

In this layout a column of the tableau is stored (mostly) contiguously in memory. Due to bitpacking, e.g., packing 64 bits into a single UInt64, the memory layout is not perfectly contiguous, but it is still optimal given that some bitwrangling is required to extract a given bit.

See also: fastrow

source
QuantumClifford.fastrowFunction

Convert a tableau to a memory layout that is fast for row operations.

In this layout a Pauli string (a row of the tableau) is stored contiguously in memory.

See also: fastrow

source
QuantumClifford.generate!Method

Generate a Pauli operator by using operators from a given the Stabilizer.

It assumes the stabilizer is already canonicalized. It modifies the Pauli operator in place, generating it in reverse, up to a phase. That phase is left in the modified operator, which should be the identity up to a phase. Returns the new operator and the list of indices denoting the elements of stabilizer that were used for the generation.

julia> ghz = S"XXXX
               ZZII
               IZZI
               IIZZ";


julia> canonicalize!(ghz)
+ XXXX
+ Z__Z
+ _Z_Z
+ __ZZ

julia> generate!(P"-ZIZI", ghz)
(- ____, [2, 4])

When the Pauli operator can not be generated by the given tableau, nothing is returned.

julia> generate!(P"XII",canonicalize!(S"ZII")) === nothing
true

julia> generate!(P"XII",canonicalize!(S"XII")) === nothing
false
source
QuantumClifford.ghzFunction

Prepare a GHZ state of n qubits.

julia> ghz()
+ XXX
+ ZZ_
+ _ZZ

julia> ghz(2)
+ XX
+ ZZ

julia> ghz(4)
+ XXXX
+ ZZ__
+ _ZZ_
+ __ZZ
source
QuantumClifford.graph_gateMethod

A helper function converting the gate indices from graphstate into a Clifford operator.

julia> s = S" XXX
              YZ_
             -_ZZ";


julia> graph, h_idx, ip_idx, z_idx = graphstate(s);


julia> gate = graph_gate(h_idx, ip_idx, z_idx, nqubits(s));


julia> apply!(s, gate) # This is now a graph state (notice you need to multiply row 1 by row 2)
+ YYZ
+ XZ_
+ _ZX

julia> canonicalize!(s) == canonicalize!(Stabilizer(graph))
true

See also: graph_gatesequence

source
QuantumClifford.graph_gatesequenceMethod

A helper function converting the gate indices from graphstate into a sequence of gates.

julia> s = S" XXX
              YZ_
             -_ZZ";


julia> graph, h_idx, ip_idx, z_idx = graphstate(s);


julia> gates = graph_gatesequence(h_idx, ip_idx, z_idx);


julia> for gate in vcat(gates...) apply!(s, gate) end


julia> s # This is now a graph state (notice you need to multiply row 1 by row 2)
+ YYZ
+ XZ_
+ _ZX

julia> canonicalize!(s) == canonicalize!(Stabilizer(graph))
true

See also: graph_gatesequence

source
QuantumClifford.graphstateMethod

Convert any stabilizer state to a graph state

Graph states are a special type of entangled stabilizer states that can be represented by a graph. For a graph $G=(V,E)$ the corresponding stabilizers are $S_v = X_v \prod_{u ∈ N(v)} Z_u$. Notice that such tableau rows contain only a single X operator. There is a set of single qubit gates that converts any stabilizer state to a graph state.

This function returns the graph state corresponding to a stabilizer and the gates that might be necessary to convert the stabilizer into a state representable as a graph.

For a tableau stab you can convert it with:

graph, hadamard_idx, iphase_idx, flips_idx = graphstate()

where graph is the graph representation of stab, and the rest specifies the single-qubit gates converting stab to graph: hadamard_idx are the qubits that require a Hadamard gate (mapping X ↔ Z), iphase_idx are (different) qubits that require an inverse Phase gate (Y → X), and flips_idx are the qubits that require a phase flip (Pauli Z gate), after the previous two sets of gates.

julia> using Graphs

julia> s = S" XXX
              ZZ_
             -_ZZ";


julia> g, h_idx, ip_idx, z_idx = graphstate(s);


julia> collect(edges(g))
2-element Vector{Graphs.SimpleGraphs.SimpleEdge{Int64}}:
 Edge 1 => 2
 Edge 1 => 3

julia> h_idx
2-element Vector{Int64}:
 2
 3

julia> ip_idx
Int64[]

julia> z_idx
1-element Vector{Int64}:
 3

The Graphs.jl library provides many graph-theory tools and the MakieGraphs.jl library provides plotting utilities for graphs.

You can directly call the graph constructor on a stabilizer, if you just want the graph and do not care about the Clifford operation necessary to convert an arbitrary state to a state representable as a graph:

julia> collect(edges( Graph(bell()) ))
1-element Vector{Graphs.SimpleGraphs.SimpleEdge{Int64}}:
 Edge 1 => 2

For a version that does not copy the stabilizer, but rather performs transformations in-place, use graphstate!. It would perform canonicalize_gott! on its argument as it finds a way to convert it to a graph state.

source
QuantumClifford.groupifyMethod

Return the full stabilizer group represented by the input generating set (a Stabilizer).

The returned object is exponentially long.

julia> groupify(S"XZ ZX")
+ __
+ XZ
+ ZX
+ YY
source
QuantumClifford.logdotMethod

Logarithm of the inner product between to Stabilizer states.

If the result is nothing, the dot inner product is zero. Otherwise the inner product is 2^(-logdot/2).

The actual inner product can be computed with LinearAlgebra.dot.

Based on (Garcia et al., 2012).

source
QuantumClifford.minimal_generating_setMethod

For a not-necessarily-minimal generating set, return the minimal generating set.

The input has to have only real phases.

julia> minimal_generating_set(S"__ XZ ZX YY")
+ XZ
+ ZX
source
QuantumClifford.normalizerMethod

Return all Pauli operators with the same number of qubits as the given Tableau t that commute with all operators in t.

julia> normalizer(T"X")
+ _
+ X
source
QuantumClifford.pauligroupMethod

Return the full Pauli group of a given length. Phases are ignored by default, but can be included by setting phases=true.

julia> pauligroup(1)
+ _
+ X
+ Z
+ Y

julia> pauligroup(1, phases=true)
+ _
+ X
+ Z
+ Y
- _
- X
- Z
- Y
+i_
+iX
+iZ
+iY
-i_
-iX
-iZ
-iY
source
QuantumClifford.pfmeasurementsMethod
pfmeasurements(frame::PauliFrame) -> Any

Returns the measurement results for each frame in the PauliFrame instance.

Relative measurements

The return measurements are relative to the reference measurements, i.e. they only say whether the reference measurements have been flipped in the given frame.

source
QuantumClifford.pfmeasurementsMethod
pfmeasurements(register::Register, frame::PauliFrame) -> Any

Takes the references measurements from the given Register and applies the flips as prescribed by the PauliFrame relative measurements. The result is the actual (non-relative) measurement results for each frame.

source
QuantumClifford.pftrajectoriesMethod
pftrajectories(
    circuit;
    trajectories,
    threads
) -> PauliFrame{Stabilizer{QuantumClifford.Tableau{Vector{UInt8}, LinearAlgebra.Adjoint{UInt64, Matrix{UInt64}}}}, Matrix{Bool}}

The main method for running Pauli frame simulations of circuits. See the other methods for lower level access.

Multithreading is enabled by default, but can be disabled by setting threads=false. Do not forget to launch Julia with multiple threads enabled, e.g. julia -t4, if you want to use multithreading.

See also: mctrajectories, petrajectories

source
QuantumClifford.pftrajectoriesMethod
pftrajectories(
    register::Register,
    circuit;
    trajectories
) -> Tuple{Register, PauliFrame{Stabilizer{QuantumClifford.Tableau{Vector{UInt8}, LinearAlgebra.Adjoint{UInt64, Matrix{UInt64}}}}, Matrix{Bool}}}

For a given Register and circuit, simulates the reference circuit acting on the register and then also simulate numerous PauliFrame trajectories. Returns the register and the PauliFrame instance.

Use pfmeasurements to get the measurement results.

source
QuantumClifford.phasesMethod

The phases of a given tableau. It is a view, i.e. if you modify this array, the original tableau caries these changes.

source
QuantumClifford.prodphaseMethod

Get the phase of the product of two Pauli operators.

Phase is encoded as F(4) in the low qubits of an UInt8.

julia> P"ZZZ"*P"XXX"
-iYYY

julia> prodphase(P"ZZZ", P"XXX")
0x03

julia> prodphase(P"XXX", P"ZZZ")
0x01
source
QuantumClifford.random_brickwork_clifford_circuitMethod

Random brickwork Clifford circuit.

The connectivity of the random circuit is brickwork in some dimensions. Each gate in the circuit is a random 2-qubit Clifford gate.

The brickwork is defined as follows: The qubits are arranged as a lattice, and lattice_size contains side length in each dimension. For example, a chain of length five will have lattice_size = (5,), and a 5×5 lattice will have lattice_size = (5, 5).

In multi-dimensional cases, gate layers act alternatively along each direction. The nearest two layers along the same direction are offset by one qubit, forming a so-called brickwork. The boundary condition is chosen as open.

source
QuantumClifford.random_pauliFunction

A random Pauli operator on n qubits.

Use nophase=false to randomize the phase. Use realphase=false to get operators with phases including ±i.

Optionally, a "flip" probability p can be provided specified, in which case each bit is set to I with probability 1-p and to X or Y or Z with probability p. Useful for simulating unbiased Pauli noise.

See also random_pauli!

source
QuantumClifford.stabilizerplotFunction

A Makie.jl recipe for pictorial representation of a tableau.

Requires a Makie.jl backend to be loaded, e.g. using CairoMakie.

Alternatively, you can use the Plots.jl plotting ecosystem, e.g. using Plots; plot(S"XXX ZZZ").

Consult the documentation for more details on visualization options.

source
QuantumClifford.stabilizerplot_axisFunction

A Makie.jl recipe for pictorial representation of a tableau.

Requires a Makie.jl backend to be loaded, e.g. using CairoMakie.

Alternatively, you can use the Plots.jl plotting ecosystem, e.g. using Plots; plot(S"XXX ZZZ").

Consult the documentation for more details on visualization options.

source
QuantumClifford.xbitMethod

Extract as a new bit array the X part of the UInt array of packed qubits of a given Pauli operator.

source
QuantumClifford.zbitMethod

Extract as a new bit array the Z part of the UInt array of packed qubits of a given Pauli operator.

source
QuantumInterface.apply!Function

In QuantumClifford the apply! function is used to apply any quantum operation to a stabilizer state, including unitary Clifford operations, Pauli measurements, and noise. Thus, this function may result in a random/stochastic result (e.g. with measurements or noise).

source
QuantumInterface.embedMethod

Embed a Pauli operator in a larger Pauli operator.

julia> embed(5, 3, P"-Y")
- __Y__

julia> embed(5, (3,5), P"-YX")
- __Y_X
source
QuantumInterface.entanglement_entropyFunction

Get bipartite entanglement entropy of a subsystem

Defined as entropy of the reduced density matrix.

It can be calculated with multiple different algorithms, the most performant one depending on the particular case.

Currently implemented are the :clip (clipped gauge), :graph (graph state), and :rref (Gaussian elimination) algorithms. Benchmark your particular case to choose the best one.

source
QuantumInterface.entanglement_entropyMethod

Get bipartite entanglement entropy by first converting the state to a graph and computing the rank of the adjacency matrix.

Based on "Entanglement in graph states and its applications".

source
QuantumInterface.expectMethod
expect(p::PauliOperator, st::AbstractStabilizer)

Compute the expectation value of a Pauli operator p on a stabilizer state st. This function will allocate a temporary copy of the stabilizer state st.

source
QuantumInterface.project!Method
project!(
    state,
    pauli::PauliOperator;
    keep_result,
    phases
) -> Tuple{MixedStabilizer, Any, Any}

Project the state of a Stabilizer on the two eigenspaces of a Pauli operator.

Assumes the input is a valid stabilizer. The projection is done inplace on that stabilizer and it does not modify the projection operator.

It returns

  • a stabilizer that might not be in canonical form
  • the index of the row where the non-commuting operator was (that row is now equal to pauli; its phase is not updated and for a faithful measurement simulation it needs to be randomized by the user)
  • and the result of the projection if there was no non-commuting operator (nothing otherwise)

If keep_result==false that result of the projection in case of anticommutation is not computed, sparing a canonicalization operation. This canonicalization operation is the only one potentially of cubic complexity. The rest of the calculations are of quadratic complexity.

If you need to measure a single qubit instead of a multiqubit Pauli operator, the faster projectX!, projectY!, and projectZ! are available.

For less boilerplate and automatic randomization of the phase use projectrand!.

Here is an example of a projection destroying entanglement:

julia> ghz = S"XXXX
               ZZII
               IZZI
               IIZZ";


julia> canonicalize!(ghz)
+ XXXX
+ Z__Z
+ _Z_Z
+ __ZZ

julia> state, anticom_index, result = project!(ghz, P"ZIII");


julia> state
+ Z___
+ Z__Z
+ _Z_Z
+ __ZZ

julia> canonicalize!(state)
+ Z___
+ _Z__
+ __Z_
+ ___Z

julia> anticom_index, result
(1, nothing)

And an example of projection consistent with the stabilizer state.

julia> s = S"ZII
             IXI
             IIY";


julia> canonicalize!(s)
+ _X_
+ __Y
+ Z__

julia> state, anticom_index, result = project!(s, P"-ZII");


julia> state
+ _X_
+ __Y
+ Z__

julia> anticom_index, result
(0, 0x02)

While not the best choice, Stabilizer can be used for mixed states, simply by providing an incomplete tableau. In that case it is possible to attempt to project on an operator that can not be generated by the provided stabilizer operators. In that case we have anticom_index==rank and result===nothing, where rank is the the new rank of the tableau, one more than the number of rows in the initial tableau. However, if keep_result was set to false, then anticom_index would stay at zero.

julia> s = S"XZI
             IZI";


julia> project!(s, P"IIX")[1]
+ X__
+ _Z_

If we had used MixedStabilizer we would have added the projector to the list of stabilizers.

julia> s = one(MixedStabilizer, 2, 3)
+ Z__
+ _Z_

julia> project!(s, P"IIX")[1]
+ Z__
+ _Z_
+ __X

However, MixedDestabilizer would be an even better choice as it has $\mathcal{O}(n^2)$ complexity instead of the $\mathcal{O}(n^3)$ complexity of *Stabilizer.

julia> s = one(MixedDestabilizer, 2, 3)
𝒟ℯ𝓈𝓉𝒶𝒷
+ X__
+ _X_
𝒳ₗ━━━
+ __X
𝒮𝓉𝒶𝒷━
+ Z__
+ _Z_
𝒵ₗ━━━
+ __Z

julia> project!(s, P"IIX")[1]
𝒟ℯ𝓈𝓉𝒶𝒷
+ X__
+ _X_
+ __Z
𝒮𝓉𝒶𝒷━
+ Z__
+ _Z_
+ __X

See the "Datastructure Choice" section in the documentation for more details.

See also: projectX!, projectY!, projectZ!, projectrand!

source
QuantumInterface.project!Method
project!(
    state::MixedStabilizer,
    pauli::PauliOperator;
    phases
) -> Tuple{MixedStabilizer, Any, Any}

When using project! on MixedStabilizer it automates some of the extra steps we encounter when implicitly using the Stabilizer datastructure to represent mixed states. Namely, it helps when the projector is not among the list of stabilizers:

julia> s = S"XZI
             IZI";


julia> ms = MixedStabilizer(s)
+ X__
+ _Z_

julia> project!(ms, P"IIY")[1]
+ X__
+ _Z_
+ __Y

Similarly to project! on Stabilizer, this function has cubic complexity when the Pauli operator commutes with all rows of the tableau. Most of the time it is better to simply use MixedDestabilizer representation.

Unlike other project! methods, this one does not allow for keep_result=false, as the correct rank or anticommutation index can not be calculated without the expensive (cubic) canonicalization operation required by keep_result=true.

See the "Datastructure Choice" section in the documentation for more details.

See also: projectX!, projectY!, projectZ!.

source
QuantumInterface.reset_qubits!Method
reset_qubits!(
    s::Stabilizer,
    newstate,
    qubits;
    phases
) -> Union{PauliOperator, Stabilizer}

Reset a given set of qubits to be in the state newstate. These qubits are traced out first, which could lead to "nonlocal" changes in the tableau.

source
QuantumInterface.tensorFunction

Tensor product between operators or tableaux.

Tensor product between CiffordOperators:

julia> tensor(CliffordOperator(sCNOT), CliffordOperator(sCNOT))
X₁ ⟼ + XX__
X₂ ⟼ + _X__
X₃ ⟼ + __XX
X₄ ⟼ + ___X
Z₁ ⟼ + Z___
Z₂ ⟼ + ZZ__
Z₃ ⟼ + __Z_
Z₄ ⟼ + __ZZ

Tensor product between PauliOperators:

julia> tensor(P"-IXYZ", P"iIXYZ")
-i_XYZ_XYZ

Tensor product between Tableaux:

julia> s = S"-XX
             +ZZ";

julia> tensor(s, s, s)
- XX____
+ ZZ____
- __XX__
+ __ZZ__
- ____XX
+ ____ZZ

julia> s = S"+XZI
             -IZI";

julia> tensor(s, s)
+ XZ____
- _Z____
+ ___XZ_
- ____Z_

See also tensor_pow.

source
QuantumInterface.tensor_powMethod

Repeated tensor product of an operators or a tableau.

For CliffordOperator:

julia> tensor_pow(CliffordOperator(sHadamard), 3)
X₁ ⟼ + Z__
X₂ ⟼ + _Z_
X₃ ⟼ + __Z
Z₁ ⟼ + X__
Z₂ ⟼ + _X_
Z₃ ⟼ + __X

For PauliOperator:

julia> tensor_pow(P"IXYZ", 2)
+ _XYZ_XYZ

For Tableaux:

julia> tensor_pow(S"Z", 4)
+ Z___
+ _Z__
+ __Z_
+ ___Z

julia> s = S"+XZI
             +IZI";

julia> tensor_pow(s, 3)
+ XZ_______
+ _Z_______
+ ___XZ____
+ ____Z____
+ ______XZ_
+ _______Z_

See also tensor.

source
QuantumInterface.traceout!Method
traceout!(
    s::Union{MixedDestabilizer, MixedStabilizer},
    qubits;
    phases,
    rank
) -> Union{Tuple{Union{MixedDestabilizer, MixedStabilizer}, Any}, MixedDestabilizer, MixedStabilizer}
source

Private API

Private Implementation Details

These functions are used internally by the library and might be drastically modified or deleted without warning or deprecation.

QuantumClifford.TableauType

Internal Tableau type for storing a list of Pauli operators in a compact form. No special semantic meaning is attached to this type, it is just a convenient way to store a list of Pauli operators. E.g. it is not used to represent a stabilizer state, or a stabilizer group, or a Clifford circuit.

source
Base.hcatMethod

Horizontally concatenates tableaux.

julia> hcat(ghz(2), ghz(2))
+ XXXX
+ ZZZZ

See also: vcat

source
Base.invMethod
inv(
    c::CliffordOperator;
    phases
) -> CliffordOperator{QuantumClifford.Tableau{Vector{UInt8}, Matrix{UInt64}}}

Inverse of a CliffordOperator

julia> inv(CliffordOperator(sCNOT))
X₁ ⟼ + XX
X₂ ⟼ + _X
Z₁ ⟼ + Z_
Z₂ ⟼ + ZZ

julia> inv(CliffordOperator(sCNOT(2, 1), 2))
X₁ ⟼ + X_
X₂ ⟼ + XX
Z₁ ⟼ + ZZ
Z₂ ⟼ + _Z

julia> inv(CliffordOperator(tHadamard))
X₁ ⟼ + Z
Z₁ ⟼ + X
source
Base.vcatMethod

Vertically concatenates tableaux.

julia> vcat(ghz(2), ghz(2))
+ XX
+ ZZ
+ XX
+ ZZ

See also: hcat

source
QuantumClifford._remove_rowcol!Method

Unexported low-level function that removes a row (by shifting all rows up as necessary)

Because MixedDestabilizer is not mutable we return a new MixedDestabilizer with the same (modified) xzs array.

Used on its own, this function will break invariants. Meant to be used with projectremove!.

source
QuantumClifford._rowmove!Method

Unexported low-level function that moves row i to row j.

Used on its own, this function will break invariants. Meant to be used in _remove_rowcol!.

source
QuantumClifford.applynoise_branchesFunction

Compute all possible new states after the application of the given noise model. Reports the probability of each one of them. Deterministic (as it reports all branches of potentially random processes), part of the Noise interface.

source
QuantumClifford.initZ!Method
initZ!(frame::PauliFrame) -> PauliFrame

Inject random Z errors over all frames and qubits for the supplied PauliFrame with probability 0.5.

Calling this after initialization is essential for simulating any non-deterministic circuit. It is done automatically by most PauliFrame constructors.

source
QuantumClifford.make_sumtype_methodFunction

``` julia> makesumtypemethod([sCNOT], :apply!, (:s,)) quote function QuantumClifford.apply!(s, g::CompactifiedGate) @cases g begin sCNOT(q1, q2) => apply!(s, sCNOT(q1, q2)) end end end

source
QuantumClifford.projectremoverand!Method

Unexported low-level function that projects a qubit and returns the result while making the tableau smaller by a qubit.

Because MixedDestabilizer is not mutable we return a new MixedDestabilizer with the same (modified) xzs array.

source
QuantumClifford.remove_column!Method

Unexported low-level function that removes a column (by shifting all columns to the right of the target by one step to the left)

Because Tableau is not mutable we return a new Tableau with the same (modified) xzs array.

source
QuantumClifford.rowdecomposeMethod

Decompose a Pauli $P$ in terms of stabilizer and destabilizer rows from a given tableaux.

For given tableaux of rows destabilizer rows $\{d_i\}$ and stabilizer rows $\{s_i\}$, there are boolean vectors $b$ and $c$ such that $P = i^p \prod_i d_i^{b_i} \prod_i s_i^{c_i}$.

This function returns p, b, c.

julia> s = MixedDestabilizer(ghz(2))
𝒟ℯ𝓈𝓉𝒶𝒷
+ Z_
+ _X
𝒮𝓉𝒶𝒷
+ XX
+ ZZ

julia> phase, destab_rows, stab_rows = QuantumClifford.rowdecompose(P"XY", s)
(3, Bool[1, 0], Bool[1, 1])

julia> im^3 * P"Z_" * P"XX" * P"ZZ"
+ XY
source
QuantumClifford.to_cpuFunction

copies the memory content of the object to CPU

You can only use this function if CUDA.jl is imported

For more advanced users to_cpu(data, element_type) will reinterpret elements of data and converts them to element_type. For example based on your CPU architecture, if working with matrices of UInt32 is faster than UInt64, you can use to_cpu(data, UInt32)

julia> using QuantumClifford: to_cpu, to_gpu

julia> using CUDA # without this import, to_cpu, to_gpu are just function

julia> stab = S"- X_Z\n+ _ZZ\n+ __Z"
- X_Z
+ _ZZ
+ __Z

julia> stab_gpu = to_gpu(stab);

julia> apply!(stab_gpu, sHadamard(1));

julia> stab_result_cpu = to_cpu(stab_gpu)
- Z_Z
+ _ZZ
+ __Z
julia> using QuantumClifford: to_cpu, to_gpu

julia> using CUDA # without this import, to_cpu, to_gpu are just function

julia> pf_gpu = to_gpu(PauliFrame(1000, 2, 2));
julia> circuit = [sMZ(1, 1), sHadamard(2), sMZ(2, 2)];
julia> pftrajectories(pf_gpu, circuit);
julia> measurements = to_cpu(pf_gpu.measurements);

See also: to_gpu

source
QuantumClifford.to_gpuFunction

copies the memory content of the object to GPU

You can only use this function if CUDA.jl is imported

For more advanced users to_gpu(data, element_type) will reinterpret elements of data and converts them to element_type. For example based on your GPU architecture, if working with matrices of UInt64 is faster than UInt32, you can use to_gpu(data, UInt64)

julia> using QuantumClifford: to_cpu, to_gpu

julia> using CUDA # without this import, to_cpu, to_gpu are just function

julia> stab = S"- X_Z\n+ _ZZ\n+ __Z"
- X_Z
+ _ZZ
+ __Z

julia> stab_gpu = to_gpu(stab);

julia> apply!(stab_gpu, sHadamard(1));

julia> stab_result_cpu = to_cpu(stab_gpu)
- Z_Z
+ _ZZ
+ __Z
julia> using QuantumClifford: to_cpu, to_gpu

julia> using CUDA # without this import, to_cpu, to_gpu are just function

julia> pf_gpu = to_gpu(PauliFrame(1000, 2, 2));
julia> circuit = [sMZ(1, 1), sHadamard(2), sMZ(2, 2)];
julia> pftrajectories(pf_gpu, circuit);
julia> measurements = to_cpu(pf_gpu.measurements);

See also: to_cpu

source