Quantum Circuits

Built-in circuits

PastaQ.ghzFunction
ghz(n::Int)

Generate a list of gates for the GHZ state.

$|\psi\rangle = (|0,0,\dots,0\rangle + |1,1,\dots,1\rangle)/\sqrt{2}$

source
PastaQ.qftFunction
qft(n::Int; inverse::Bool = false)

Generate a list of gates for the quantum fourier transform circuit on $n$ qubits.

source

Assembling circuits

PastaQ.gatelayerFunction
gatelayer(gatename::AbstractString, n::Int; kwargs...)

Create a uniform layer containing n identical quantum gates, idenfitied by gatename. If additional parameteres are provided, they are identically added to all gates.

gatelayer("H",3)
# 3-element Vector{Tuple{String, Int64}}:
#  ("H", 1)
#  ("H", 2)
#  ("H", 3)

gatelayer("X",1:2:5)
# 3-element Vector{Tuple{String, Int64}}:
#  ("X", 1)
#  ("X", 3)
#  ("X", 5)
source
gatelayer(gatename::AbstractString, bonds::Vector{Vector{Int}}; kwargs...)

Create a uniform layer of multi-qubit gates over a set of bonds.

gatelayer("CX", [(j,j+1), j=1:2:5])
# 3-element Vector{Tuple{String, Tuple{Int64, Int64}}}:
#  ("CX", (1, 2))
#  ("CX", (3, 4))
#  ("CX", (5, 6))
source
PastaQ.randomlayerFunction
randomlayer(
  gatename::AbstractString,
  support::Union{Int, Vector{<:Tuple}, AbstractRange};
  rng = Random.GLOBAL_RNG,
  kwargs...
)

Generate a random layer built out of a set one or two qubit gates If support::Int = n, generates $n$ single-qubit gates gatename. If support::Vector=bonds, generates a set of two-qubit gates on the couplings contained in support.

randomlayer("Ry", 3)
# 3-element Vector{Any}:
#  ("Ry", 1, (θ = 0.5029516521736083,))
#  ("Ry", 2, (θ = 2.5324545964433693,))
#  ("Ry", 3, (θ = 2.0510824561219523,))
source
randomlayer(
  gatenames::Vector{<:AbstractString},
  support::Union{Vector{<:Int}, AbstractRange, Vector{<:Tuple}};
  rng=Random.GLOBAL_RNG,
  weights::Union{Nothing,Vector{Float64}} = ones(length(gatenames)) / length(gatenames),
  kwargs...,
)

Generate a random layer built out of one or two qubit gates, where gatenames is a set of possible gates to choose from. By default, each single gate is sampled uniformaly over this set. If weights are provided, each gate is sampled accordingly.

randomlayer(["X","Y","Z"], 3)
# 3-element Vector{Any}:
#  ("Y", 1)
#  ("Y", 2)
#  ("X", 3)
source

Random circuits

PastaQ.randomcircuitFunction
randomcircuit(
  coupling_sequence::Vector;
  depth::Int = 1,
  twoqubitgates::Union{String,Vector{String}}="RandomUnitary",
  onequbitgates::Union{Nothing,String,Vector{String}}=nothing,
  layered::Bool=true,
  rng=Random.GLOBAL_RNG)

Build a circuit with given depth, where each layer consists of a set of two-qubit gates applied on pairs of qubits in according to a set of coupling_sequences. Each layer also contains $n$ single-qubit gates. IN both cases, the chosen gates are passed as keyword arguments onequbitgates and twoqubitgates.

The default configurations consists of two-qubit random Haar unitaries, and no single-qubit gates.

If layered = true, the object returned in a Vector of circuit layers, rather than the full collection of quantum gates.

source
randomcircuit(n::Int; kwargs...)

One-dimensional random quantum circuit:

randomcircuit(4; depth = 2, twoqubitgates = "CX", onequbitgates = "Ry")
# [("CX", (1, 2)),
#  ("CX", (3, 4)),
#  ("Ry", 1, (θ = 0.52446,)),
#  ("Ry", 2, (θ = 3.01059,)),
#  ("Ry", 3, (θ = 0.25144,)),
#  ("Ry", 4, (θ = 1.93356,))]
# [("CX", (2, 3)),
#  ("Ry", 1, (θ = 2.15460,)),
#  ("Ry", 2, (θ = 2.52480,)),
#  ("Ry", 3, (θ = 1.85756,)),
#  ("Ry", 4, (θ = 0.02405,))]
source
randomcircuit(size::Tuple; rotated::Bool = false, kwargs...)

Two-dimensional random quantum circuit on a square lattice. If rotated = true, use rotated lattice of 45 degrees.

source

Executing quantum circuits

PastaQ.buildcircuitFunction
buildcircuit(
  hilbert::Vector{<:Index},
  circuit::Union{Tuple, Vector{<:Any}};
  noise::Union{Nothing, Tuple, NamedTuple} = nothing
)
buildcircuit(M::Union{MPS,MPO,ITensor}, args...; kwargs...)

Compile a circuit from a lazy representation into a vector of ITensor. For example, a gate element of circuit, ("gn", (i,j)) is turned into a rrank-4 tensor corresponding to the i and j element of hilbert.

If noise is passed, the corresponding Kraus operator are inserted appropriately after the gates in the circuit.

source
PastaQ.runcircuitFunction
runcircuit(circuit::Any; kwargs...)

Execute quantum circuit (see below).

source
runcircuit(hilbert::Vector{<:Index}, circuit::Tuple; kwargs...)

Execute quantum circuit on Hilbert space hilbert (see below).

source
runcircuit(hilbert::Vector{<:Index}, circuit::Vector;
           full_representation::Bool = false,
           process::Bool = false,
           noise = nothing,
           kwargs...)

runcircuit(M::Union{MPS, MPO, ITensor}, circuit::Union{Tuple, AbstractVector};
           full_representation::Bool = false, noise = nothing, kwargs...)

Run the circuit corresponding to a list of quantum gates on a system of $n$ qubits, with input Hilbert space hilbert. The specific method of this general function is specified by the keyword arguments process and noise.

By default (process = false and noise = nothing), runcircuit returns an MPS wavefunction corresponding to the contraction of each quantum gate in circuit with the zero product state

\[|\psi\rangle = U_M\dots U_2 U_1|0,0,\dots,0\rangle\]

If process = true, the output is the MPO corresponding to the full unitary circuit:

\[U = U_M\dots U_2 U_1\]

If noise is set to a given input noise model (in lazy representation), Kraus operators are added to each gate in the circuit, and the output is the MPO density operator given by the contraction of the noisy circuit with input zero state:

\[\rho = \mathcal{E}(|0,\dots,0\rangle\langle0,\dots,0|)\]

Finally, if both noise = ... and proces = true, the output is the full quantum channel, which we by default represent with its Choi matrix:

\[\Lambda = (1+\mathcal{E})|\Phi\rangle\langle\Phi|^{\otimes n}\]

If full_representation = true, the contraction is performed without approximation, leading to an output object whose size scales exponentiall with $n$

source
runcircuit(
  M::Union{MPS, MPO, ITensor},
  circuit::Vector{<:Vector{<:ITensor}};
  (observer!)=nothing,
  move_sites_back_before_measurements::Bool=false,
  noise = nothing,
  outputlevel = 1,
  outputpath = nothing,
  savestate = false,
  print_metrics = [],
  kwargs...)

Apply a quantum circuit to an input state M, where the circuit is built out of a sequence of layers of quantum gates. The input state may be an MPS wavefunction $|\psi\rangle$, an MPO density operator $ρ$ (or unitary operator $U$), etc.

By feeding a "layered" circuit, we can enable measurement and keep track of metrics as a function of the circuit's depth.

Other than the keyword arguments of the high-level interface, here we can provide:

  • (observer!): observer object (from Observers.jl).
  • outputlevel = 1: amount of printing during calculation.
  • outputpath = nothing: if set, save observer on file.
  • savestate = false: if true, save the MPS/MPO on file.
  • print_metrics = []: the metrics in the observed to print at each depth.`
source
runcircuit(
  M::Union{MPS,MPO},
  circuit_tensors::Vector{<:ITensor};
  apply_dag=nothing,
  cutoff=1e-15,
  maxdim=10_000,
  svd_alg="divide_and_conquer",
  move_sites_back::Bool=true,
  kwargs...)

Apply a set of "gate" tensors (alredy in the form of ITensor) to an input state M, with options:

  • apply_dag = nothing: whether to perform conjugate evolution.
  • cutoff = 1e-15: truncation cutoff in SVD.
  • maxdim = 10_000: maximum bond dimension at SVD trunctions.
  • svd_alg = "divide_and_conquer": SVD algorithm (see ITensors.jl).
  • move_sites_back = true: move sites back after long-range gate.

By default, apply_dag = nothing and the interface is dictated by the input state, and whether or not the vector of ITensor containins rank-3 noisy tensors (i.e. Kraus operators).

For an input MPS $|\psi_0\rangle$, with a unitary circuit, the output is

\[|\psi\rangle = U_M\dots U_2 U_1|\psi_0\rangle\]

while for noisy circuits:

\[\rho = \mathcal{E}(|\psi_0\rangle\langle\psi_0|)\]

For an input MPO $\rho_0$, the output is

\[\rho = U_M\dots U_2 U_1 \rho_0 U^\dagger_1 U^\dagger_2,\dots,U^\dagger_M\]

for unitary circuits, and $\rho = \mathcal{E}(\rho_0)$ for noisy circuits.

source
PastaQ.choimatrixFunction
choimatrix(circuit::Vector{<:Any}; kwargs...) =
choimatrix(sites::Vector{<:Index}, circuit::Vector{<:Any}; kwargs...)
choimatrix(sites::Vector{<:Index}, circuit_tensors::Vector{<:ITensor};
           full_representation = false,kwargs...)

Compute the Choi matrix for a noisy channel

\[\Lambda = (1+\mathcal{E})|\Phi\rangle\langle\Phi|^{\otimes n}\]

If full_representation = true, the contraction is performed without approximation, leading to an output object whose size scales exponentiall with $n$

source