RandomMeas
Documentation for RandomMeas.jl: The randomized measurement toolbox in Julia
- Data acquisition
- Postprocessing (excluding classical shadows)
- Postprocessing for classical shadows
- Simulating quantum circuits
- Additional useful routines for ITensor
- Examples
Data acquisition
Data acquisition types
RandomMeas.AbstractMeasurementSetting
— TypeAbstractMeasurementSetting
An abstract type representing a general measurement setting. Concrete implementations (e.g. LocalUnitaryMeasurementSetting
, ComputationalBasisMeasurementSetting
) should subtype this.
RandomMeas.ComputationalBasisMeasurementSetting
— TypeComputationalBasisMeasurementSetting
A struct representing computational basis measurement settings for quantum systems. This setting uses the computational basis, so that each local unitary is by construction simply the identity operator.
Fields
N::Int
: Number of sites (qubits).local_unitary::Vector{ITensor}
: A vector of N identity ITensors.site_indices::Vector{Index{Int64}}
: A vector of site indices (length N).
Constraints
N == length(site_indices)
.
RandomMeas.ComputationalBasisMeasurementSetting
— MethodComputationalBasisMeasurementSetting(ms::ComputationalBasisMeasurementSetting;
N=ms.N,
local_unitary=ms.local_unitary,
site_indices=ms.site_indices)
Make a new ComputationalBasisMeasurementSetting
by copying fields from ms
, but overriding any that you pass by keyword.
RandomMeas.LocalMeasurementSetting
— TypeLocalMeasurementSetting
An abstract type for measurement settings that correspond to local (i.e. single qubit/site) measurements.
RandomMeas.LocalUnitaryMeasurementSetting
— TypeLocalUnitaryMeasurementSetting(N, local_unitary, site_indices)
A measurement setting where each qubit is specified by a single-qubit rotation. Rotates from the computational basis into the measurement basis.
Fields
N::Int
: Number of sites (qubits).local_unitary::Vector{ITensor}
: A vector ofN
ITensors representing the local unitary basis rotations.site_indices::Vector{Index{Int64}}
: A vector of site indices of lengthN
.
Constraints
N == length(local_unitary) == length(site_indices)
.- Each ITensor in
local_unitary
has exactly two indices:- One unprimed (
site_indices[i]
) - One primed (
prime(site_indices[i])
).
- One unprimed (
RandomMeas.LocalUnitaryMeasurementSetting
— MethodLocalUnitaryMeasurementSetting(ms::LocalUnitaryMeasurementSetting;
N=ms.N,
local_unitary=ms.local_unitary,
site_indices=ms.site_indices)
Make a new LocalUnitaryMeasurementSetting
by copying fields from ms
, but overriding any that you pass by keyword.
RandomMeas.MeasurementData
— TypeMeasurementData{T}
A container for measurement data and settings obtained in actual or simulated quantum experiments.
Fields
N::Int
: Number of sites (qubits).NM::Int
: Number of measurements per setting.measurement_results::Array{Int, 2}
: A 2D array of binary measurement results with dimensions(NM, N)
.measurement_setting::T
: A measurement setting of typeT
(subtype ofAbstractMeasurementSetting
) ornothing
.
Type Parameter
T
: The type of the measurement setting, constrained toUnion{Nothing, AbstractMeasurementSetting}
.
Usage
Typically constructed via the provided constructors.
RandomMeas.MeasurementGroup
— TypeMeasurementGroup{T}
A container for a group of measurement data objects used in actual or simulated quantum experiments.
Fields
N::Int
: Number of sites (qubits).NU::Int
: Number of measurement data objects.NM::Int
: Number of measurements per setting.measurements::Vector{MeasurementData{T}}
: A vector of measurement data objects.
Type Parameter
T
: The type of the measurement setting for each measurement data object, constrained toUnion{Nothing, AbstractMeasurementSetting}
.
Usage
Typically constructed via one of the provided constructors.
Example
# Assume setting1 and setting2 are valid measurement settings
data1 = MeasurementData(results1; measurement_setting=setting1)
data2 = MeasurementData(results2; measurement_setting=setting2)
group = MeasurementGroup([data1, data2])
RandomMeas.MeasurementProbability
— TypeMeasurementProbability{T}
A container for measurement probabilities and settings obtained either estimated from measurement data or directly computed from quantum states.
Fields
N::Int
: Number of sites (qubits).measurement_probability::ITensor
: An ITensor representing Born probability.measurement_setting::T
: A measurement setting of typeT
ornothing
.site_indices::Vector{Index{Int64}}
: A vector of site indices (length N).
Type Parameter
T
: The type of measurement setting, constrained toUnion{Nothing, AbstractMeasurementSetting}
.
Usage
Constructed either from measurement data or directly from quantum states.
RandomMeas.ShallowUnitaryMeasurementSetting
— TypeShallowUnitaryMeasurementSetting
A struct representing measurement settings which is, for each qubit, specified through a single qubit rotation, rotating from the computational basis into the measurement basis.
Fields
N::Int
: Number of sites (qubits).- 'K::Int`: Number of gates that creates the shallow_unitary
localunitary::Vector{ITensor}
: A vector of Ngates representing the shallow unitarysite_indices::Vector{Index{Int64}}
: A vector of site indices of length N.
Constructor
Creates a ShallowUnitaryMeasurementSetting
object after validating that:
- The length of
local_unitary
equalsK
- The length of
site_indices
equalsN
.
RandomMeas.ShallowUnitaryMeasurementSetting
— MethodShallowUnitaryMeasurementSetting(ms::ShallowUnitaryMeasurementSetting;
N=ms.N,
K=ms.K,
local_unitary=ms.local_unitary,
site_indices=ms.site_indices)
Make a new ShallowUnitaryMeasurementSetting
by copying fields from ms
, but overriding any that you pass by keyword.
Data acquisition routines
RandomMeas.ComputationalBasisMeasurementSetting
— MethodComputationalBasisMeasurementSetting(N; site_indices=nothing)
Create a ComputationalBasisMeasurementSetting
for N
sites. This setting corresponds to measurement in the computational basis.
Arguments:
N::Int
: Number of sites (qubits).site_indices::Union{Vector{Index{Int64}}, Nothing}
(optional): Site indices. Ifnothing
, they are automatically generated.
Returns:
- A
ComputationalBasisMeasurementSetting
object.
Example:
setting = ComputationalBasisMeasurementSetting(4)
RandomMeas.LocalUnitaryMeasurementSetting
— MethodLocalUnitaryMeasurementSetting(local_unitary_array; site_indices=nothing)
Create a LocalUnitaryMeasurementSetting
object from an N × 2 × 2
array of unitary matrices.
Arguments:
local_unitary_array::Array{ComplexF64, 3}
: AnN × 2 × 2
array of unitary matrices.site_indices::Union{Vector{Index{Int64}}, Nothing}
(optional): Site indices. Ifnothing
, they are automatically generated.
Returns:
- A
LocalUnitaryMeasurementSetting
object.
Example:
unitary_array = rand(ComplexF64, 4, 2, 2)
setting = LocalUnitaryMeasurementSetting(unitary_array)
RandomMeas.LocalUnitaryMeasurementSetting
— MethodLocalUnitaryMeasurementSetting(N; site_indices=nothing, ensemble="Haar")
Create a LocalUnitaryMeasurementSetting
object by randomly sampling local unitary operators.
Arguments:
N::Int
: Number of sites (qubits).site_indices::Union{Vector{Index}, Nothing}
(optional): Site indices. Ifnothing
, they are automatically generated.ensemble::String
: Type of random unitary ("Haar"
,"Pauli"
,"Identity"
).
Returns:
- A
LocalUnitaryMeasurementSetting
object.
Example:
setting = LocalUnitaryMeasurementSetting(4, ensemble="Haar")
RandomMeas.ShallowUnitaryMeasurementSetting
— MethodShallowUnitaryMeasurementSetting(N, depth; site_indices=nothing)
Create a ShallowUnitaryMeasurementSetting
object by generating a random quantum circuit.
Arguments:
N::Int
: Number of sites (qubits).depth::Int
: Depth of the random circuit.site_indices::Union{Vector{Index{Int64}}, Nothing}
(optional): Site indices. Ifnothing
, they are automatically generated.
Returns:
- A
ShallowUnitaryMeasurementSetting
object.
Example:
setting = ShallowUnitaryMeasurementSetting(4, 3)
RandomMeas.export_LocalUnitaryMeasurementSetting
— MethodExport the unitary in a LocalUnitaryMeasurementSetting object to an .npz file with a single field: local_unitary.
Arguments:
ms::LocalUnitaryMeasurementSetting
: The measurement settings to export.filepath::String
: Path to the output .npz file.
RandomMeas.get_rotation
— Functionget_rotation(site_index, ensemble="Haar")
Generate a single-qubit unitary sampled from a specified ensemble.
Arguments:
site_index::Index{Int64}
: Site index.ensemble::String
: Type of unitary ensemble ("Haar"
,"Pauli"
,"Identity"
).
Returns:
- An
ITensor
representing the unitary transformation.
Example:
U = get_rotation(site_index, "Pauli")
RandomMeas.import_LocalUnitaryMeasurementSetting
— MethodImport unitary from an .npz file and create a LocalUnitaryMeasurementSetting object.
Arguments:
filepath::String
: Path to the input .npz file.site_indices::Union{Vector{Index{Int64}}, Nothing}
: Optional site indices. If not provided, they will be generated.
Returns:
- A LocalUnitaryMeasurementSettings object.
RandomMeas.reduce_to_subsystem
— Methodreduce_to_subsystem(settings, subsystem)
Reduce a LocalMeasurementSetting
object to a specified subsystem.
Arguments:
settings::LocalMeasurementSetting
: The original measurement settings object.subsystem::Vector{Int}
: A vector of site indices (1-based) specifying the subsystem to retain.
Returns:
- A new
LocalMeasurementSetting
object corresponding to the specified subsystem.
Example:
reduced_setting = reduce_to_subsystem(full_setting, [1, 3])
RandomMeas.MeasurementData
— MethodMeasurementData(ψ::Union{MPO, MPS}, NM::Int; mode::String = "MPS/MPO", measurement_setting::Union{LocalUnitaryMeasurementSetting, ComputationalBasisMeasurementSetting, ShallowUnitaryMeasurementSetting} = nothing)
Returns a MeasurementData
object by sampling NM
projective measurements from the quantum state ψ
.
Arguments
ψ::Union{MPO, MPS}
: The quantum state represented as a Matrix Product Operator (MPO) or Matrix Product State (MPS).NM::Int
: The number of measurement shots to simulate for each setting.mode::String
(optional): Specifies the simulation method. Options:"dense"
: Uses the dense representation."MPS/MPO"
(default): Uses tensor network methods for memory efficiency.
measurement_setting
(optional): A measurement setting object (if not provided, defaults to computational basis measurements).
Returns
A MeasurementData
object with the corresponding measurement results and setting.
RandomMeas.MeasurementData
— MethodMeasurementData(measurement_results::Array{Int, 2}; measurement_setting::Union{T, Nothing} = nothing)
Creates a MeasurementData
object by inferring the dimensions of the measurement results and validating the provided setting.
Arguments
measurement_results::Array{Int, 2}
: A 2D array of binary measurement results with shape(NM, N)
.measurement_setting::Union{T <: AbstractMeasurementSetting, Nothing}
(optional): Measurement setting ornothing
if not provided.
Returns
A MeasurementData
object with inferred dimensions and validated setting.
Throws
AssertionError
: If the dimensions ofmeasurement_results
andmeasurement_setting
are inconsistent.
Examples
# With measurement setting
setting = LocalUnitaryMeasurementSetting(4, ensemble="Haar")
results = rand(1:2, 10, 4)
data_with_setting = MeasurementData(results; measurement_setting=setting)
# Without measurement setting
data_without_setting = MeasurementData(rand(1:2, 10, 4))
RandomMeas.MeasurementData
— MethodMeasurementData(measurement_probability::MeasurementProbability{T}, NM::Int) where T <: Union{Nothing, AbstractMeasurementSetting}
Returns a MeasurementData
object by sampling NM
projective measurements based on the provided measurement probability.
Arguments
measurement_probability::MeasurementProbability
: A container with the measurement probability (an ITensor) and associated settings.NM::Int
: The number of projective measurements to sample.
Returns
A MeasurementData
object with dimensions inferred from the measurement probability.
RandomMeas.export_MeasurementData
— Methodexport_measurement_data(data::MeasurementData, filepath::String)
Exports measurement data to a .npz
file.
Arguments
data::MeasurementData
: The measurement data object containing measurement results and optionally a LocalUnitaryMeasurementSetting setting.filepath::String
: The file path where the data will be exported.
Details
- The
measurement_results
are exported directly as they are. - If
measurement_setting
is provided, the associatedlocal_unitaries
are extracted, reshaped, and included in the export.
Notes
- The exported
.npz
file will contain:"measurement_results"
: A 2D array of shape(NM, N)
, where:NM
: Number of measurements per setting.N
: Number of qubits/sites.
"local_unitaries"
(optional): A 4D array of shape(N, 2, 2)
representing the unitary transformations for each site.
Example
# Create MeasurementData object
data = MeasurementData(measurement_results; measurement_setting=measurement_setting)
# Export to a file
export_measurement_data(data, "exported_data.npz")
RandomMeas.import_MeasurementData
— Methodimport_measurement_data(filepath::String; predefined_setting=nothing, site_indices=nothing)
Imports measurement results and optional measurement settings from an archive file.
Arguments
filepath::String
: Path to the.npz
file containing the measurement results and optionally local unitaries. The file should contain at least a fieldmeasurement_results
(2D binary array of shape(NM, N)
), and optionally a fieldlocal_unitaries
(local unitaries as a Nx2x2 array).predefined_setting
(optional): A predefinedMeasurementSetting
object. If provided, this will be used instead of the file's setting.site_indices
(optional): A vector of site indices to be used when constructingLocalUnitaryMeasurementSetting
from the fieldlocal_unitaries
(only relevant if predefined_setting is not provided). If not provided, the default site indices will be generated internally.
Returns
A MeasurementData
object containing the imported results and settings.
Examples
# Import with predefined settings
setting = LocalUnitaryMeasurementSetting(local_unitaries; site_indices=siteinds("Qubit", 5))
data_with_setting = import_measurement_data("data.npz"; predefined_setting=setting)
# Import with site indices provided
data_with_indices = import_measurement_data("data.npz"; site_indices=siteinds("Qubit", 5))
# Import without any additional options
data = import_measurement_data("data.npz")
RandomMeas.reduce_to_subsystem
— Methodreduce_to_subsystem(data::MeasurementData{T}, subsystem::Vector{Int}) where T <: Union{Nothing, LocalMeasurementSetting}
Reduce a MeasurementData
object to a specified subsystem, preserving the measurement setting type if available.
Arguments
data::MeasurementData{T}
: The original measurement data object, whereT
is eithernothing
or a subtype ofLocalMeasurementSetting
.subsystem::Vector{Int}
: A vector of site indices (1-based) specifying the subsystem to retain. Each index must be between 1 anddata.N
.
Returns
A new MeasurementData{T}
object with:
- The measurement results reduced from dimensions
(NM, N)
to(NM, |subsystem|)
. - The measurement setting reduced accordingly (if one is provided), or remaining as
nothing
.
Example
# Suppose `data` is a MeasurementData object with N = 4.
# To retain only sites 1 and 3:
reduced_data = reduce_to_subsystem(data, [1, 3])
RandomMeas.MeasurementGroup
— MethodMeasurementGroup(ψ::Union{MPO, MPS}, NU::Int, NM::Int, depth::Int; mode::String = “MPS/MPO”, progress_bar::Bool=false)
::MeasurementGroup{ShallowUnitaryMeasurementSetting}
Construct a MeasurementGroup from a quantum state ψ
by generating NU
shallow measurement settings and simulating NM
measurements per unitary.
Arguments
ψ::Union{MPO, MPS}
: The quantum state.NU::Int
: Number of measurement data objects to generate.NM::Int
: Number of measurements per setting.depth::Int
: Circuit depth for shallow settings.mode::String
: Simulation mode; defaults to “MPS/MPO”.progress_bar
::Bool: Whether to show a progress bar.
Returns
A MeasurementGroup{ShallowUnitaryMeasurementSetting} object.
RandomMeas.MeasurementGroup
— MethodMeasurementGroup(ψ::Union{MPO, MPS}, NU::Int, NM::Int; mode::String = “MPS/MPO”, progress_bar::Bool=false)
::MeasurementGroup{LocalUnitaryMeasurementSetting}
Construct a MeasurementGroup from a quantum state ψ
by generating NU
local measurement settings and simulating NM
projective measurements per setting.
Arguments
ψ::Union{MPO, MPS}
: The quantum state.NU::Int
: Number of measurement data objects to generate.NM::Int
: Number of measurements per setting.mode::String
: Simulation mode; defaults to “MPS/MPO”.progress_bar::Bool
: Whether to show a progress bar.
Returns
A MeasurementGroup{LocalUnitaryMeasurementSetting} object.
RandomMeas.MeasurementGroup
— MethodMeasurementGroup(measurements::Vector{MeasurementData{T}}) where {T <: Union{Nothing, AbstractMeasurementSetting}}
Construct a MeasurementGroup
object by inferring dimensions from a vector of MeasurementData
objects.
Arguments
measurements::Vector{MeasurementData{T}}
: A vector ofMeasurementData
objects.
Returns
A MeasurementGroup
object with:
N
: Inferred from the first element (assumed consistent across all elements).NU
: Number of measurement data objects.NM
: Inferred from the first element.measurements
: The provided vector.
Example
setting1 = LocalUnitaryMeasurementSetting(4, ensemble="Haar")
results1 = rand(1:2, 10, 4)
data1 = MeasurementData(results1; measurement_setting=setting1)
setting2 = LocalUnitaryMeasurementSetting(4, ensemble="Haar")
results2 = rand(1:2, 10, 4)
data2 = MeasurementData(results2; measurement_setting=setting2)
group = MeasurementGroup([data1, data2])
RandomMeas.MeasurementGroup
— MethodMeasurementGroup(ψ::Union{MPO, MPS}, measurement_settings::Vector{AbstractMeasurementSetting}, NM::Int; mode::String = “MPS/MPO”, progress_bar::Bool=false)
::MeasurementGroup{T} where T <: AbstractMeasurementSetting
Construct a MeasurementGroup from a quantum state ψ
by generating NU
local measurement settings and simulating NM
projective measurements per setting.
Arguments
ψ::Union{MPO, MPS}
: The quantum state.measurement_settings::Vector{AbstractMeasurementSetting}
: A vector with measurement settingsNM::Int
: Number of measurements per setting.mode::String
: Simulation mode; defaults to “MPS/MPO”.progress_bar::Bool
: Whether to show a progress bar.
Returns
A MeasurementGroup{T} object.
RandomMeas.export_MeasurementGroup
— Methodexport_MeasurementGroup(group::MeasurementGroup{T}, filepath::String)
Export a MeasurementGroup object to an NPZ file.
Arguments
group::MeasurementGroup{T}
: A MeasurementGroup object where each MeasurementData may have its own measurement setting of type T (with T <: Union{Nothing, LocalUnitaryMeasurementSetting, ComputationalBasisMeasurementSetting}).filepath::String
: The file path where the NPZ file will be written.
Details
- The measurement results from each MeasurementData object (each of shape (NM, N)) are stacked into a 3D array of shape (NU, NM, N), where NU is the number of MeasurementData objects.
- The measurement settings are exported as a Complex array of size NU x N x 2 x 2 if present.
RandomMeas.import_MeasurementGroup
— Methodimport_MeasurementGroup(filepath::String; predefined_settings=nothing, site_indices=nothing) -> MeasurementGroup
Import a MeasurementGroup object from an NPZ file.
Arguments
filepath::String
: The path to the NPZ file containing the exported MeasurementGroup data.predefined_settings
(optional): A vector of predefined measurement settings (one per MeasurementData object). If provided, its length must equal the exported NU.site_indices
(optional): A vector of N site indices to use when reconstructing the measurement setting. If not provided, default site indices are generated usingsiteinds("Qubit", N)
.
Returns
A MeasurementGroup object with:
- Measurement results reconstructed from a 3D array of shape (NU, NM, N).
- A measurement setting for each MeasurementData object reconstructed from a 4D array of shape (NU, N, 2, 2) if present, or taken from
predefined_settings
if provided.
RandomMeas.reduce_to_subsystem
— Methodreduce_to_subsystem(
group::MeasurementGroup{T},
subsystem::Vector{Int}
)::MeasurementGroup{T} where T <: LocalMeasurementSetting
Reduce a MeasurementGroup object (with
LocalUnitaryMeasurementSetting`) to a specified subsystem.
Arguments
group::MeasurementGroup{LocalUnitaryMeasurementSetting}
: The original measurement data object.subsystem::Vector{Int}
: A vector of site indices (1-based) specifying the subsystem to retain.
Returns
A new MeasurementGroup
object corresponding to the specified subsystem.
Postprocessing (excluding classical shadows)
RandomMeas.get_XEB
— Methodget_XEB(ψ::MPS, measurement_data::MeasurementData)
Return the linear cross-entropy for the measurement results in measurement_data
, with respect to a theory state ψ
.
Arguments:
ψ::MPS
: The theoretical state to compare against.measurement_data::MeasurementData
: The measurement data object containing results and settings.
Returns:
The linear cross-entropy as a Float64
.
RandomMeas.get_fidelity
— Functionget_fidelity(
group_1::MeasurementGroup,
group_2::MeasurementGroup,
subsystem::Vector{Int} = collect(1:group_1.N)
)
Compute the fidelity of two quantum states Tr(ρ1 ρ2)/SROOT(Tr(ρ1^2),Tr(ρ2^2)) from measurement data by averaging the overlap of measurement results.
Arguments
group_1::MeasurementGroup
: Measurement data for the first state.group_2::MeasurementGroup
: Measurement data for the second state.subsystem::Vector{Int}
(optional): A vector of site indices specifying the subsystem to retain. Defaults to the full system.
Returns
- The computed fidelity.
RandomMeas.get_h_tensor
— Functionget_h_tensor(s::Index, s_prime::Index) -> ITensor
Construct the Hamming tensor for given indices.
Arguments
s::Index
: Unprimed site index.s_prime::Index
: Primed site index.
Returns
Hamming_tensor::ITensor
: The Hamming tensor connectings
ands_prime
.
Method
- Initializes an
ITensor
with indicess
ands_prime
. - Assigns values to represent the Hamming distance operation:
- Diagonal elements are set to
1.0
. - Off-diagonal elements are set to
-0.5
.
- Diagonal elements are set to
RandomMeas.get_overlap
— Functionget_overlap(
group_1::MeasurementGroup,
group_2::MeasurementGroup,
subsystem::Vector{Int} = collect(1:group_1.N);
apply_bias_correction::Bool = false
)
Compute the overlap of two quantum states from measurement data by averaging the overlap of measurement results.
Arguments
group_1::MeasurementGroup
: Measurement data for the first state.group_2::MeasurementGroup
: Measurement data for the second state.subsystem::Vector{Int}
(optional): A vector of site indices specifying the subsystem to retain. Defaults to the full system.apply_bias_correction::Bool
(optional): Whether to apply bias correction for the overlap. Defaults tofalse
.
Returns
- The computed overlap (or purity if
group_1 == group_2
and bias correction is applied).
RandomMeas.get_overlap
— Methodget_overlap(
data_1::MeasurementData,
data_2::MeasurementData
)
Compute the overlap between two quantum states for a single measurement setting.
Arguments
data_1::MeasurementData
: Measurement Data for the first state, with dimensions(NM, N)
.data_2::MeasurementData
: Measurement Data for the second state, with dimensions(NM, N)
.
Returns
- The computed overlap for the single measurement setting.
RandomMeas.get_overlap
— Methodget_overlap(prob1::MeasurementProbability, prob2::MeasurementProbability) -> Float64
Compute the weighted overlap \2^N sum_s (-2)^{-D[s,s']}P(s)P(s')]
by sequentially applying the Hamming tensor to each qubit index and contracting with the second probability tensor.
Arguments
prob1::MeasurementProbability
: The first Born probability tensor representing quantum staterho1
.prob2::MeasurementProbability
: The second Born probability tensor representing quantum staterho2
.
Returns
weighted_overlap::Float64
: The computed traceTr(rho1 rho2)
scaled appropriately..
Example
using ITensors
# Assume prob1 and prob2 are predefined MeasurementProbabilities
overlap = get_overlap(prob1, prob2)
println("Overlap: ", overlap)
RandomMeas.get_purity
— Functionget_purity(group::Measurementgroup, subsystem::Vector{Int} = collect(1:group.N))
Compute the purity of a quantum state from measurement data by averaging the overlap of measurement results.
Arguments
group::MeasurementGroup
: Measurement data containing the results and settings of randomized measurements.subsystem::Vector{Int}
(optional): A vector of site indices specifying the subsystem to retain. Defaults to the full system.
Returns
- The computed purity for the specified subsystem.
Postprocessing for classical shadows
Classical shadow types
RandomMeas.AbstractShadow
— TypeAbstractShadow
An abstract type representing a general classical shadow. Concrete subtypes should implement specific shadow methodologies, such as factorized or dense shadows.
RandomMeas.DenseShadow
— TypeDenseShadow
A struct representing a dense classical shadow (a 2^N x 2^N matrix), stored as a single ITensor with 2N indices.
Fields
shadow_data::ITensor
: An ITensor with 2N indices representing the dense shadow.N::Int
: Number of sites (qubits).site_indices::Vector{Index{Int64}}
: A vector of site indices (length N).
Constructor
DenseShadow(shadow_data::ITensor, N::Int, site_indices::Vector{Index{Int64}})
validates that:
site_indices
has length N.shadow_data
has exactly 2N indices.- The set of unprimed indices in
shadow_data
matchessite_indices
. - The set of primed indices in
shadow_data
matchesmap(prime, site_indices)
.
RandomMeas.FactorizedShadow
— TypeFactorizedShadow
A struct representing a factorized classical shadow which can be represented as a tensor product of single qubit shadows.
Fields
shadow_data::Vector{ITensor}
: A vector of ITensors (each 2×2) representing the shadow for each qubit/site.N::Int
: Number of qubits/sites.site_indices::Vector{Index{Int64}}
: A vector of site indices (length N).
Constructor
FactorizedShadow(shadow_data::Vector{ITensor}, N::Int, site_indices::Vector{Index{Int64}})
validates that:
- The length of
shadow_data
andsite_indices
equalsN
. - Each ITensor in
shadow_data
has exactly two indices, which include the corresponding unprimed and primed site index.
RandomMeas.ShallowShadow
— TypeShallowShadow
A struct representing a shallow classical shadow, stored as a MPO ITensor object.
Fields
shadow_data::MPOr
: An MPO representing the shallow shadow.N::Int
: Number of sites (qubits).site_indices::Vector{Index{Int64}}
: A vector of site indices (length N).
Constructor
ShallowShadow(shadow_data::MPO, N::Int, site_indices::Vector{Index{Int64}})
validates that:
site_indices
has length N.shadow_data
has exactly N Tensors, and the site indices match with site_indices
Classical shadow routines
RandomMeas.get_expect_shadow
— Methodget_expect_shadow(O::MPO, shadows::AbstractArray{<:AbstractShadow}; compute_sem::Bool = false)
Compute the average expectation value of an MPO operator O
using an array of shadow objects.
Arguments
O::MPO
: The MPO operator whose expectation value is to be computed.shadows::AbstractArray{<:AbstractShadow}
: An array of shadow objects (of any shape) over which the expectation values are computed.compute_sem::Bool
(optional): Iftrue
, also compute the standard error of the mean (SEM). Default isfalse
.
Returns
- If
compute_sem
isfalse
, returns the average expectation value. - If
compute_sem
istrue
, returns a tuple(mean, sem)
, wheremean
is the average expectation value andsem
is the standard error.
Example
mean_val = get_expect_shadow(O, shadows)
mean_val, sem_val = get_expect_shadow(O, shadows; compute_sem=true)
RandomMeas.get_expect_shadow
— Methodget_expect_shadow(O::MPO, shadow::AbstractShadow)
Compute the expectation value of an MPO operator O
using a single shadow object.
Arguments:
O::MPO
: The MPO operator for which the expectation value is computed.shadow::AbstractShadow
: A shadow object, either dense, factorized, or shallow.
Returns
The expectation value as a scalar.
Example
val = get_expect_shadow(O, shadow)
RandomMeas.get_trace_moment
— Methodget_trace_moment(shadows::Array{<:AbstractShadow, 2}, kth_moment::Int; O::Union{Nothing, MPO}=nothing)
Compute a single trace moment from an array of AbstractShadow
objects.
Arguments
shadows::Array{<:AbstractShadow, 2}
: An array of shadow objects with dimensions(n_ru, n_m)
, wheren_ru
is the number of random unitaries andn_m
is the number of measurements.kth_moment::Int
: The momentk
to compute (e.g.,k = 1, 2, ...
).O::Union{Nothing, MPO}
(optional): If provided, computesTr[O * ρ^k]
; otherwise, computesTr[ρ^k]
(default:nothing
).
Returns
The computed trace moment as a scalar.
Example
moment1 = get_trace_moment(shadows, 1)
moment2 = get_trace_moment(shadows, 2; O=my_operator)
RandomMeas.get_trace_moment
— Methodget_trace_moment(shadows::Vector{<:AbstractShadow}, kth_moment::Int; O::Union{Nothing, MPO}=nothing)
Wrapper function. Compute a single trace moment for a vector of shadow objects by reshaping the vector into a 2D array.
Arguments
shadows::Vector{<:AbstractShadow}
: A vector of shadow objects.kth_moment::Int
: The moment orderk
to compute.O::Union{Nothing, MPO}
(optional): An MPO observable.
Returns
The computed trace moment as a scalar.
Example
moment = get_trace_moment(shadows_vector, 2; O=my_operator)
RandomMeas.get_trace_moments
— Methodget_trace_moments(shadows::Array{<:AbstractShadow, 2}, kth_moments::Vector{Int}; O::Union{Nothing, MPO}=nothing)
Wrapper function. Compute multiple trace moments from an array of shadow objects.
Arguments
shadows::Array{<:AbstractShadow, 2}
: An array of shadow objects with dimensions(n_ru, n_m)
.kth_moments::Vector{Int}
: A vector of moment orders.O::Union{Nothing, MPO}
(optional): An MPO observable; if provided, computesTr[O * ρ^k]
for each moment (default:nothing
).
Returns
A vector of trace moments corresponding to each moment in kth_moments
.
Example
moments = get_trace_moments(shadows_array, [1, 2, 3])
RandomMeas.get_trace_moments
— Methodget_trace_moments(shadows::Vector{<:AbstractShadow}, kth_moments::Vector{Int}; O::Union{Nothing, MPO}=nothing)
Wrapper function. Compute multiple trace moments from a vector of shadow objects by reshaping the vector into a 2D array.
Arguments
shadows::Vector{<:AbstractShadow}
: A vector of shadow objects.kth_moments::Vector{Int}
: A vector of moment orders.O::Union{Nothing, MPO}
(optional): An MPO observable.
Returns
A vector of trace moments.
Example
moments = get_trace_moments(shadows_vector, [1, 2, 3])
RandomMeas.get_trace_product
— Methodget_trace_product(shadows::AbstractShadow...; O::Union{Nothing, MPO}=nothing)
Compute the product of multiple shadow objects and return its trace or expectation value.
If O
is nothing
, returns the trace of the product: trace(shadow₁ * shadow₂ * ... * shadowₙ). If O
is provided, returns the expectation value computed by: getexpectshadow(O, shadow₁ * shadow₂ * ... * shadowₙ).
Arguments
shadows...
: A variable number of shadow objects.O::Union{Nothing, MPO}
(optional): An MPO observable.
Returns
The trace of the product if O
is nothing
, or the expectation value if O
is provided.
Example
result = get_trace_product(shadow1, shadow2, shadow3)
result_with_O = get_trace_product(shadow1, shadow2, shadow3; O=my_operator)
RandomMeas.multiply
— Methodmultiply(shadow1::AbstractShadow, shadow2::AbstractShadow)
Multiply two shadow objects of the same concrete type.
Arguments:
shadow1::AbstractShadow
: The first shadow object.shadow2::AbstractShadow
: The second shadow object.
Returns
A new shadow object representing the product.
Throws
An ArgumentError
if the types of shadow1
and shadow2
do not match.
Example
prod_shadow = multiply(shadow1, shadow2)
RandomMeas.partial_trace
— Methodpartial_trace(shadows::AbstractArray{<:AbstractShadow}, subsystem::Vector{Int})
Compute the partial trace for each shadow in a collection of shadows.
Arguments
shadows::AbstractArray{<:AbstractShadow}
: A collection of shadow objects (vector or 2D array).subsystem::Vector{Int}
: A vector of site indices (1-based) specifying the subsystem to retain.
Returns
An array of shadows reduced to the specified subsystem, with the same dimensions as the input array.
RandomMeas.partial_trace
— Methodpartial_trace(shadow::AbstractShadow, subsystem::Vector{Int}; assume_unit_trace::Bool=false)
Compute the partial trace of a shadow object over the complement of the specified subsystem.
Arguments
shadow::AbstractShadow
: The shadow object.subsystem::Vector{Int}
: A vector of site indices (1-based) specifying the subsystem to retain.assume_unit_trace::Bool
(optional): Iftrue
, assumes the shadow has unit trace (default:false
). This can speed up the calculation for factorized shadows (as the trace of "traced out" qubits is not computed)/
Returns
A new shadow object reduced to the specified subsystem.
Example
reduced_shadow = partial_trace(shadow, [1, 3])
RandomMeas.partial_transpose
— Methodpartial_transpose(shadows::AbstractArray{<:AbstractShadow}, subsystem::Vector{Int})
Compute the partial transpose for each shadow in a collection.
Arguments
shadows::AbstractArray{<:AbstractShadow}
: A collection (vector or 2D array) of shadow objects.subsystem::Vector{Int}
: A vector of site indices (1-based) specifying the subsystem(s) over which the transpose is to be performed.
Returns
An array of shadow objects with the partial transpose applied, preserving the input dimensions.
RandomMeas.partial_transpose
— Methodpartial_transpose(shadow::AbstractShadow, subsystem::Vector{Int})
Compute the partial transpose of a shadow object over the specified subsystem(s). This operation is analogous to QuTiP's partial transpose method.
Arguments
shadow::AbstractShadow
: The shadow object for which the partial transpose is computed.subsystem::Vector{Int}
: A vector of site indices (1-based) specifying the subsystem(s) over which to perform the transpose.
Returns
A new shadow object that is the partial transpose of the input.
Example
transposed_shadow = partial_transpose(shadow, [2, 4])
RandomMeas.trace
— Methodtrace(shadows::AbstractArray{<:AbstractShadow})
Compute the trace for each shadow in a collection of shadow objects.
Arguments:
shadows::AbstractArray{<:AbstractShadow}
: A collection (vector, matrix, etc.) of shadow objects.
Returns
An array of scalar trace values corresponding to each shadow, with the same dimensions as the input.
RandomMeas.trace
— Methodtrace(shadow::AbstractShadow)
Compute the trace of a shadow object.
Arguments:
shadow::AbstractShadow
: A shadow object.
Returns
The trace of the shadow object as a scalar.
Example
t = trace(shadow)
RandomMeas.FactorizedShadow
— MethodFactorizedShadow(measurement_results::Vector{Int}, local_unitary::Vector{ITensor};
G::Vector{Float64} = fill(1.0, length(local_unitary)))
Construct a FactorizedShadow
object from raw measurement results and unitary transformations.
Arguments
measurement_results::Vector{Int}
: Vector of binary measurement results for each qubit/site.local_unitary::Vector{ITensor}
: Vector of local unitary transformations applied during the measurement.- G::Vector{Float64}
(optional): Vector of
G` values for measurement error correction (default: 1.0 for all sites).
Returns
A FactorizedShadow
object.
RandomMeas.convert_to_dense_shadow
— Methodconvert_to_dense_shadow(factorized_shadow::FactorizedShadow)
Convert a FactorizedShadow
object into a DenseShadow
object.
Arguments
factorized_shadow::FactorizedShadow
: TheFactorizedShadow
object to convert.
Returns
A DenseShadow
object with the combined ITensor.
RandomMeas.get_expect_shadow
— Methodget_expect_shadow(O::MPO, shadow::FactorizedShadow)
Compute the expectation value of an MPO operator O
using a factorized shadow.
Arguments:
O::MPO
: The MPO operator for which the expectation value is computed.shadow::FactorizedShadow
: A factorized shadow object.
Returns:
The expectation value as a ComplexF64
(or Float64
if purely real).
RandomMeas.get_factorized_shadows
— Methodget_factorized_shadows(measurement_data::MeasurementData{LocalUnitaryMeasurementSetting};
G::Vector{Float64} = fill(1.0, measurement_data.N))
Compute factorized shadows for all measurement results in the provided MeasurementData
.
Arguments
measurement_data::MeasurementData{LocalUnitaryMeasurementSetting}
: Measurement data object containing measurement results and settings.G::Vector{Float64}
(optional): Vector ofG
values for measurement error correction (default: 1.0 for all sites).
Returns
A Vector of NM FactorizedShadow
objects with dimensions.
RandomMeas.get_factorized_shadows
— Methodget_factorized_shadows(measurement_group::MeasurementGroup{LocalUnitaryMeasurementSetting};
G::Vector{Float64} = fill(1.0, measurement_group.N))
Compute factorized shadows for all measurement results in the provided MeasurementGroup
.
Arguments
measurement_group::MeasurementGroup{LocalUnitaryMeasurementSetting}
: Measurement data object containing measurement results and settings.G::Vector{Float64}
(optional): Vector ofG
values for measurement error correction (default: 1.0 for all sites).
Returns
A Array of NU*NM FactorizedShadow
objects with dimensions.
RandomMeas.multiply
— Methodmultiply(shadow1::FactorizedShadow, shadow2::FactorizedShadow)
Multiply two FactorizedShadow
objects element-wise.
Arguments
shadow1::FactorizedShadow
: The firstFactorizedShadow
object.shadow2::FactorizedShadow
: The secondFactorizedShadow
object.
Returns
A new FactorizedShadow
object representing the element-wise product of the two inputs.
Notes
- Both
shadow1
andshadow2
must have the same number of qubits/sites.
RandomMeas.partial_trace
— Methodpartial_trace(shadow::FactorizedShadow, subsystem::Vector{Int}; assume_unit_trace::Bool = false)
Compute the partial trace of a FactorizedShadow
object over the complement of the specified subsystem.
Arguments
shadow::FactorizedShadow
: The factorized shadow to compute the partial trace for.subsystem::Vector{Int}
: A vector of site indices (1-based) specifying the subsystem to retain.assume_unit_trace::Bool
(optional): Iftrue
, assumes all traced-out tensors have unit trace and skips explicit computation (default:false
).
Returns
A new FactorizedShadow
object reduced to the specified subsystem.
Notes
- If
assume_unit_trace
istrue
, avoids explicit trace computation for efficiency. - If
assume_unit_trace
isfalse
, computes the traces of all tensors outside the subsystem and multiplies their product into the remaining tensors. - Issues a warning if the trace product deviates significantly from 1 when
assume_unit_trace
isfalse
.
RandomMeas.partial_transpose
— Methodpartial_transpose(shadow::FactorizedShadow, subsystem::Vector{Int})::FactorizedShadow
Compute the partial transpose of a FactorizedShadow over the specified subsystem by swapping, for each site, the unprimed and primed indices using the swapind
function. This function returns views of the underlying ITensors, avoiding unnecessary data duplication.
Arguments
shadow::FactorizedShadow
: The factorized classical shadow.subsystem::Vector{Int}
: A vector of 1-based site indices on which to perform the partial transpose.
Returns
A new FactorizedShadow with the specified sites partially transposed.
RandomMeas.trace
— Methodtrace(shadow::FactorizedShadow)
Compute the trace of a FactorizedShadow
object.
Arguments
shadow::FactorizedShadow
: TheFactorizedShadow
object whose trace is to be computed.
Returns
The trace of the shadow as a Float64
or ComplexF64
.
Notes
- The function computes the product of the traces of individual tensors in the factorized shadow.
RandomMeas.DenseShadow
— MethodDenseShadow(measurement_data::MeasurementData{LocalUnitaryMeasurementSetting}; G::Vector{Float64} = fill(1.0, size(measurement_data.N, 2)))
Construct a DenseShadow
object from a MeasurementDataObject
Arguments
- `measurement_data::MeasurementData{LocalUnitaryMeasurementSetting}:
G::Vector{Float64}
(optional): Vector of G values to account for measurement errors (default: 1.0 for all sites).
Returns
A DenseShadow
object.
RandomMeas.DenseShadow
— MethodDenseShadow(measurement_probability::MeasurementProbability; G::Vector{Float64} = fill(1.0, length(u)))
Construct a DenseShadow
object from a precomputed probability tensor.
Arguments
P::ITensor
: Probability tensor representing measurement results.u::Vector{ITensor}
: Vector of local unitary transformations.G::Vector{Float64}
(optional): Vector of G values to account for measurement errors (default: 1.0 for all sites).
Returns
A DenseShadow
object.
RandomMeas.get_dense_shadows
— Methodget_dense_shadows(measurement_group::MeasurementGroup{LocalUnitaryMeasurementSetting};
G::Vector{Float64} = fill(1.0, N),
number_of_ru_batches::Int = NU)
Compute dense shadows for the provided measurement data in batches.
Arguments
measurement_group::MeasurementGroup{LocalUnitaryMeasurementSetting}
: Measurement data object.G::Vector{Float64}
(optional): Vector of G values for robustness (default: 1.0 for all sites).number_of_ru_batches::Int
(optional): Number of random unitary batches (default:NU
).
Returns
A Vector of DenseShadow
objects.
RandomMeas.get_expect_shadow
— Methodget_expect_shadow(O::MPO, shadow::DenseShadow)
Compute the expectation value of an MPO operator O
using a dense shadow.
Arguments:
O::MPO
: The MPO operator for which the expectation value is computed.shadow::DenseShadow
: A dense shadow object.
Returns:
The expectation value as a ComplexF64
(or Float64
if purely real).
RandomMeas.multiply
— Methodmultiply(shadow1::DenseShadow, shadow2::DenseShadow)
Compute the trace product of two dense shadows.
Arguments
shadow1::DenseShadow
: The first dense shadow.shadow2::DenseShadow
: The second dense shadow.
Returns
A new DenseShadow
object that represents the trace product of the two input shadows.
Notes
- The shadows must have the same site indices (
ξ
) and number of qubits (N
).
RandomMeas.partial_trace
— Methodpartial_trace(shadow::DenseShadow, subsystem::Vector{Int})
Compute the partial trace of a DenseShadow
object over the complement of the specified subsystem.
Arguments
shadow::DenseShadow
: The dense shadow to compute the partial trace for.subsystem::Vector{Int}
: A vector of site indices (1-based) specifying the subsystem to retain.
Returns
A new DenseShadow
object reduced to the specified subsystem.
RandomMeas.partial_transpose
— Methodpartial_transpose(shadow::DenseShadow, subsystem::Vector{Int})::DenseShadow
Compute the partial transpose of a DenseShadow over the specified subsystem by swapping, for each site, the unprimed index with its primed partner. This is done using the swapind
function, which returns a view of the underlying ITensor.
Arguments
shadow::DenseShadow
: The dense classical shadow.subsystem::Vector{Int}
: A vector of 1-based site indices on which to perform the partial transpose.
Returns
A new DenseShadow with the specified sites partially transposed.
RandomMeas.trace
— Methodtrace(shadow::DenseShadow)
Compute the trace of a DenseShadow
object.
Arguments
shadow::DenseShadow
: TheDenseShadow
object whose trace is to be computed.
Returns
The trace of the shadow as a Float64
or ComplexF64
.
Notes
- The function contracts the
ξ
andξ'
indices of the shadow's ITensor.
RandomMeas.ShallowShadow
— MethodShallowShadow(measurement_results::Vector{Int}, local_unitary::Vector{ITensor};
G::Vector{Float64} = fill(1.0, length(local_unitary)))
Construct a ShallowSShadow
object from raw measurement results and unitary transformations.
Arguments
measurement_results::Vector{Int}
: Vector of binary measurement results for each qubit/site.local_unitary::Vector{ITensor}
: Vector of local unitary transformations applied during the measurement.
Returns
A ShallowShadow
object.
RandomMeas.apply_map
— Methodapply_map(map::MPO,state::MPO,s::Vector{Index{Int64}},ξ::Vector{Index{Int64}})
Apply a map map((s,s')→(ξ,ξ')) on a state of indices (ξ,ξ')
RandomMeas.get_depolarization_map
— Methodget_depolarization_map(depolarization_mps::MPS,s::Vector{Index{Int64}},ξ::Vector{Index{Int64}})
returns a shallow map \mathcal{M} parametrization by a depolarizationmps c(\nu) where the state is depolarized over partition \A{ u} with probability c(\nu)=1
RandomMeas.get_depolarization_map
— Methodget_depolarization_map(depolarization_mps::MPS,s::Vector{Index{Int64}},ξ::Vector{Index{Int64}})
returns a shallow map \mathcal{M} parametrization by a depolarizationmps c(\nu) where the state is depolarized over partition \A{ u} with probability c(\nu)=1
RandomMeas.get_expect_shadow
— Methodget_expect_shadow(O::MPO, measurement_data::MeasurementData{ShallowUnitaryMeasurementSetting}, inverse_shallow_map::Vector{ITensor},s::Vector{Index{Int64}},ξ::Vector{Index{Int64}})
Compute the expectation value of an MPO operator O
from a shallow MeasurementData and and inverse shallow_map
Returns:
The expectation value as a ComplexF64
(or Float64
if purely real).
RandomMeas.get_expect_shadow
— Methodget_expect_shadow(O::MPO, shadow::ShallowShadow)
Compute the expectation value of an MPO operator O
using a shallow shadow.
Arguments:
O::MPO
: The MPO operator for which the expectation value is computed.shadow::ShallowShadow
: A factorized shadow object.
Returns:
The expectation value as a ComplexF64
(or Float64
if purely real).
RandomMeas.get_shallow_depolarization_mps
— Methodget_shallow_depolarization_mps(group::MeasurementGroup{ShallowUnitaryMeasurementSetting})
TBW
RandomMeas.get_shallow_shadows
— Methodget_shallow_shadows(measurement_data::MeasurementData{ShallowUnitaryMeasurementSetting}, inverse_shallow_map::Vector{ITensor},s::Vector{Index{Int64}},ξ::Vector{Index{Int64}})
Construct a ShallowSShadow
object from MeasurementData
Arguments
Returns
A ShallowShadow
object.
Simulating quantum circuits
RandomMeas.apply_depo_channel
— Methodapply_depo_channel(ρ::MPO, p::Vector{Float64})
Apply a local depolarization channel to an MPO by modifying each site tensor according to the depolarization probability.
For each site, the channel acts as:
ρ[i] → (1 - p[i]) * ρ[i] + (p[i] / 2) * (ρ[i] * δ(s, s') * δ(s, s'))
where δ(s, s') is the delta tensor that contracts the site index with its primed counterpart.
Arguments
ρ::MPO
: The input Matrix Product Operator representing a density matrix.p::Vector{Float64}
: A vector of depolarization probabilities, one per site.
Returns
An MPO with the depolarization channel applied on each site.
RandomMeas.apply_depo_channel
— Methodapply_depo_channel(ψ::MPS, p::Vector{Float64})
Apply the local depolarization channel to an MPS by converting it to an MPO density matrix (using the outer product) and then applying the depolarization channel.
Arguments
ψ::MPS
: The input Matrix Product State representing a pure state.p::Vector{Float64}
: A vector of depolarization probabilities, one per site.
Returns
An MPO representing the depolarized density matrix.
RandomMeas.random_Pauli_layer
— Methodrandom_Pauli_layer(ξ::Vector{Index{Int64}}, p::Vector{Float64})
Construct a layer of random single-qubit Pauli operations to simulate local depolarization. Upon avereraging, this corresponds to the local depolarization channel with strength p.
For each qubit (with index i), a random Pauli operation is applied with the following probabilities:
- With probability 1 - 3p_i/4: No operation is applied (the qubit remains unchanged).
- With probability p_i/4} each: Apply the X, Y, or Z gate.
Here, p_i is the depolarization probability for qubit i.
Arguments
ξ::Vector{Index{Int64}}
: A vector of ITensor indices representing the qubit sites.p::Vector{Float64}
: A vector of depolarization probabilities (one per qubit).
Returns
A vector of ITensors representing the applied Pauli gates. If no gate is applied on a site (with probability 1 - 3p_i/4, that site is omitted from the returned circuit.
Example
circuit = random_Pauli_layer(siteinds("Qubit", 5), 0.05 * ones(5))
RandomMeas.random_circuit
— Methodrandom_circuit(ξ::Vector{Index{Int64}}, depth::Int64)
Create a random circuit of the given depth. The function returns a vector of ITensors, each representing a gate in the circuit.
- If
depth == 0
, a single-qubit random unitary is applied to each site. - For
depth > 0
, the circuit is built layer-by-layer:- On odd layers, random two-qubit gates are applied on sites 1-2, 3-4, etc.
- On even layers, random two-qubit gates are applied on sites 2-3, 4-5, etc.
Arguments
ξ::Vector{Index{Int64}}
: A vector of site indices for the qubits.depth::Int64
: The depth of the circuit (non-negative integer).
Returns
A vector of ITensors representing the random circuit gates.
Example
circuit = random_circuit(siteinds("Qubit", 10), 3)
RandomMeas.random_magnetic_field_layer
— Methodrandom_magnetic_field_layer(ξ::Vector{Index{Int64}}, p::Vector{Float64})
Construct a layer of random Rz gates representing a random magnetic field along the z-axis.
For each qubit i, a random rotation Rz is applied with a rotation angle drawn uniformly from [0, 2 pi pi). This gives an average rotation angle of pi pi on each site.
Arguments
ξ::Vector{Index{Int64}}
: A vector of ITensor indices corresponding to the qubit sites.p::Vector{Float64}
: A vector of parameters (one per qubit) that set the scale of the rotation angles.
Returns
A vector of ITensors representing the random Rz gates applied to each qubit.
Example
circuit = random_magnetic_field_layer(siteinds("Qubit", 5), 0.1 * ones(5))
Additional useful routines for ITensor
RandomMeas.flatten
— Methodflatten(O::Union{MPS, MPO, Vector{ITensor}})
Flatten a Matrix Product State (MPS), Matrix Product Operator (MPO), or a vector of ITensors into a single ITensor by sequentially multiplying the constituent tensors.
Arguments
O
: An MPS, MPO, or vector of ITensors to be flattened.
Returns
An ITensor representing the product of the individual tensors in O
.
Example
A = random_mps(siteinds("Qubit", 5))
flatA = flatten(A)
RandomMeas.get_Born_MPS
— Methodget_Born_MPS(ρ::MPO)
Construct the Born probability vector as an MPS from an MPO representation of a density matrix ρ.
This function computes the Born probability vector P(s) = ⟨s|ρ|s⟩, where |s⟩ is a basis state. It does so by contracting each tensor of the MPO ρ with appropriate delta tensors that enforce equality between the unprimed and primed indices. The result is returned as an MPS that represents the Born probabilities over the computational basis.
Arguments
ρ::MPO
: A Matrix Product Operator representing the density matrix.
Returns
An MPS representing the Born probability vector.
Example
P = get_Born_MPS(ρ)
RandomMeas.get_Born_MPS
— Methodget_Born_MPS(ψ::MPS)
Construct the Born probability vector P(s) = |ψ(s)|² as an MPS from an MPS representation ψ.
This function computes the probability for each computational basis state by contracting each tensor of the MPS ψ with its complex conjugate, using appropriate delta tensors to enforce index equality. The resulting MPS represents the Born probability distribution of the state.
Arguments
ψ::MPS
: A matrix product state representing a pure quantum state.
Returns
An MPS representing the Born probability vector, where each tensor P[i] corresponds to the probability contribution at site i.
Example
P = get_Born_MPS(ψ)
RandomMeas.get_average_mps
— Methodget_average_mps(ψ_list::Vector{MPS}, χ::Int64, nsweeps::Int64)
Approximate the average state σ from a collection of MPS using a DMRG-like algorithm.
The algorithm finds an MPS ψ (with maximum bond dimension χ) that approximates the average state σ = Average(ψ_list). To monitor convergence, it tracks a cost function defined as:
cost_function = ⟨ψ|ψ⟩ - ⟨ψ|σ⟩ - ⟨σ|ψ⟩,
which is equivalent to (||σ - ψ||² - ⟨σ|σ⟩).
Arguments
ψ_list::Vector{MPS}
: A vector of MPS objects representing individual quantum states.χ::Int64
: The desired maximum bond dimension for the averaged MPS.nsweeps::Int64
: The number of sweeps (iterations) to perform in the DMRG-like algorithm.
Returns
An MPS representing the approximate average state with bond dimension χ.
Example
avg_state = get_average_mps(ψ_list, 20, 10)
RandomMeas.get_entanglement_spectrum
— Methodget_entanglement_spectrum(ψ::MPS, NA::Int64)
Compute the entanglement spectrum for the bipartition defined by the first NA
sites of the MPS ψ
.
This function first creates a copy of ψ
and orthogonalizes it up to site NA
. Then, it performs an SVD on the tensor at site NA
:
- If
NA > 1
, the SVD is taken with respect to the indices corresponding to the link between sitesNA-1
andNA
and the physical index at siteNA
. - If
NA == 1
, only the physical index at siteNA
is used.
The returned object spec
contains the singular values (which are related to the Schmidt coefficients) for the bipartition.
Arguments
ψ::MPS
: A matrix product state representing a pure quantum state.NA::Int64
: The number of sites from the left that define the subsystem for which the entanglement spectrum is computed.
Returns
spec
: An ITensor containing the singular values from the SVD at siteNA
.
Example
spectrum = get_entanglement_spectrum(ψ, 3)
RandomMeas.get_selfXEB
— Methodget_selfXEB(ψ::MPS)
Compute the self-XEB (cross-entropy benchmarking) metric for a pure state represented as an MPS.
The self-XEB is defined as:
selfXEB = 2^N * ∑ₛ |ψ(s)|⁴ - 1
where the sum is over all computational basis states s and N is the number of sites (qubits). This function first computes the Born probability MPS from ψ, then calculates the inner product of the probability MPS with itself, scales the result by 2^N, and finally subtracts 1.
Arguments
ψ::MPS
: A Matrix Product State representing a pure quantum state.
Returns
A scalar (Float64) representing the self-XEB value.
Example
x = get_selfXEB(ψ)
RandomMeas.get_siteinds
— Methodget_siteinds(ψ::Union{MPS, MPO})
Retrieve the site indices for a quantum state represented as an MPS or MPO.
Arguments
ψ
: The quantum state, which can be either a Matrix Product State (MPS) or a Matrix Product Operator (MPO).
Returns
A vector of site indices corresponding to the quantum state ψ
.
Example
ξ = get_siteinds(ψ)
RandomMeas.get_trace
— Methodget_trace(ρ::MPO)
Compute the trace of a Matrix Product Operator (MPO) ρ by contracting each tensor with a delta function that equates its unprimed and primed indices.
Arguments
ρ::MPO
: A matrix product operator representing a quantum state or operator.
Returns
A scalar representing the trace of ρ.
Example
t = get_trace(ρ)
RandomMeas.get_trace_moment
— Functionget_trace_moment(ψ::Union{MPS, MPO}, k::Int, subsystem::Vector{Int}=collect(1:length(ψ)))
Compute the kth trace moment of the reduced density matrix for a given subsystem of a quantum state.
For a pure state (MPS) and when the subsystem is contiguous starting from site 1, the function computes the entanglement spectrum of the bipartition defined by the last site in subsystem
and returns the kth moment of the squared Schmidt coefficients. Otherwise, the function reduces the state to the specified subsystem.
For k = 2 (purity), it returns the squared norm (which is equivalent to tr(ρ²)). For k > 2, it computes ρ^k via repeated application of the apply
function (with a cutoff) and returns the trace of the resulting tensor.
Arguments
ψ::Union{MPS, MPO}
: The quantum state, represented as an MPS (for pure states) or MPO (for mixed states).k::Int
: The moment order to compute (must be an integer ≥ 1).subsystem::Vector{Int}
(optional): A vector of site indices (1-based) specifying the subsystem to retain. Defaults to all sites.
Returns
A scalar (Float64) representing the kth trace moment of the reduced density matrix.
Example
moment = get_trace_moment(ψ, 3, [1, 2, 3]; partial_transpose=false)
RandomMeas.get_trace_moment
— Methodget_trace_moment(spec::ITensor, k::Int)
Compute the kth moment of the entanglement spectrum represented by the ITensor spec
.
The function assumes that spec
is a square ITensor whose diagonal elements correspond to the singular values (Schmidt coefficients) of a reduced density matrix. The kth moment is computed as:
pk = ∑ₗ (spec[l, l]^(2*k))
which effectively computes the sum over the kth powers of the squared singular values.
Arguments
spec::ITensor
: A square ITensor representing the entanglement spectrum.k::Int
: The moment order to compute (must be an integer ≥ 1).
Returns
A scalar (Float64) corresponding to the kth moment.
Example
moment = get_trace_moment(spec, 2)
RandomMeas.get_trace_moments
— Functionget_trace_moments(ψ::Union{MPS, MPO}, k_vector::Vector{Int}, subsystem::Vector{Int}=collect(1:length(ψ)))
Compute a vector of trace moments for a quantum state ψ
over a specified subsystem.
For each moment order k
in k_vector
, the function computes the kth trace moment of the reduced density matrix obtained by applying reduce_to_subsystem(ψ, subsystem)
.
Arguments
ψ::Union{MPS, MPO}
: The quantum state, represented as an MPS (for pure states) or an MPO (for mixed states).k_vector::Vector{Int}
: A vector of integer moment orders (each ≥ 1) for which the trace moments are computed.subsystem::Vector{Int}
(optional): A vector of site indices (1-based) specifying the subsystem to consider. Defaults to all sites.
Returns
A vector of scalars, each being the kth trace moment corresponding to the entries of k_vector
.
Example
moments = get_trace_moments(ψ, [1, 2, 3], [1, 2, 3])
RandomMeas.partial_transpose
— Methodpartial_transpose(ρ::MPO, subsystem::Vector{Int})
Compute the partial transpose of an MPO over the sites specified by subsystem
.
For each site index in the MPO:
- If the index is in the
subsystem
, the tensor is transposed by swapping its unprimed and primed indices usingswapind
. - Otherwise, the tensor is left unchanged (multiplied by 1.0 for type consistency).
Arguments
ρ::MPO
: A Matrix Product Operator representing a density matrix.subsystem::Vector{Int}
: A vector of site indices (1-based) over which to apply the transpose.
Returns
An MPO in which the tensors corresponding to the sites in subsystem
have been partially transposed.
Example
ρT = partial_transpose(ρ, [2, 3])
RandomMeas.reduce_to_subsystem
— Functionreduce_to_subsystem(ρ::MPO, subsystem::Vector{Int64})
Compute the reduced density matrix (as an MPO) for a specified subsystem.
Arguments
ρ::MPO
: A Matrix Product Operator representing the full density matrix.subsystem::Vector{Int64}
: A vector of site indices (1-based) specifying the subsystem to retain.
Returns
An MPO representing the reduced density matrix over the sites specified in subsystem
.
Example
ρ_reduced = reduce_to_subsystem(ρ, [2, 3])
RandomMeas.reduce_to_subsystem
— Methodreduce_to_subsystem(ψ::MPS, subsystem::Vector{Int64})
Compute the reduced density matrix for a pure state represented by the MPS ψ
over the specified subsystem.
This function first constructs the density matrix by taking the outer product of ψ
with itself, and then applies the MPO version of reduce_to_subsystem
to obtain the reduced density matrix for the sites specified in subsystem
.
Arguments
ψ::MPS
: A Matrix Product State representing a pure quantum state.subsystem::Vector{Int64}
: A vector of site indices (1-based) specifying the subsystem to retain.
Returns
An MPO representing the reduced density matrix over the specified subsystem.
Example
ρ_sub = reduce_to_subsystem(ψ, [2, 3])
Examples
Classical shadows
Quantum benchmark
Entanglement
Analyzing the experimental data of Brydges et al, Science 2019
Surface code and the measurement of the topological entanglement entropy
Mixed-state entanglement: The p3-PPT condition and batch shadows