Functions
Detailed explanations or additional notes about all functions.
Measurement functions
BlueTangle.sample
— Functionsample(state::AbstractVectorS, shots::Int) -> Vector
Sample outcomes from a quantum state vector based on the probability distribution.
state
: A sparse quantum state vector.shots
: Number of samples to be taken.
Returns a vector of sampled outcomes.
BlueTangle.get_probs_from_sample
— Functionget_probs_from_sample(sample::Vector, N::Int) -> (Vector, Vector)
Convert a sample of outcomes into probabilities.
sample
: A vector of sampled outcomes.N
: Number of qubits.
Returns a tuple of vectors: the first vector contains outcomes, and the second vector contains corresponding probabilities.
BlueTangle.expect
— FunctionAlias:
expect(m::Measurement, qubit::Int) -> Float64
expect(m::Measurement) -> Vector
Calculate the expectation value from a measurement.
m
: AMeasurement
object.qubit
: (Optional) Specific qubit for which to calculate the expectation value.
Returns the expectation value(s).
Alias for density matrix:
expect(state::AbstractVectorS, op::QuantumOps) -> Float64
expect(state::AbstractVectorS, op_str::String, qubit::Int) -> Float64
expect(state::AbstractVectorS, matrix::sa.SparseMatrixCSC) -> Float64
expect(state::AbstractVectorS, op_str::String) -> Vector{Float64}
Alias for state vector:
expect(rho::sa.SparseMatrixCSC, op::QuantumOps) -> Float64
expect(rho::sa.SparseMatrixCSC, op_str::String, qubit::Int) -> Float64
expect(rho::sa.SparseMatrixCSC, matrix::sa.SparseMatrixCSC) -> Float64
expect(rho::sa.SparseMatrixCSC, op_str::String) -> Vector{Float64}
Calculate the expectation value for quantum states or density matrices given an operator. This function has several forms depending on the input parameters:
expect(state::AbstractVectorS, op::QuantumOps)
: Computes the expectation value for a quantum state vector with a given operator.expect(rho::sa.SparseMatrixCSC, op::QuantumOps)
: Computes the expectation value for a density matrix with a given operator.expect(state::AbstractVectorS, op_str::String, qubit::Int)
: Computes the expectation value for a specific qubit in a quantum state vector with an operator specified as a string.expect(rho::sa.SparseMatrixCSC, op_str::String, qubit::Int)
: Computes the expectation value for a specific qubit in a density matrix with an operator specified as a string.expect(state::AbstractVectorS, op_str::String)
: Computes the expectation values for all qubits in a quantum state vector given an operator as a string.expect(rho::sa.SparseMatrixCSC, op_str::String)
: Computes the expectation values for all qubits in a density matrix given an operator as a string.expect(state::AbstractVectorS, matrix::sa.SparseMatrixCSC)
: Computes the expectation value using a sparse matrix representation of an operator for a state vector.expect(rho::sa.SparseMatrixCSC, matrix::sa.SparseMatrixCSC)
: Computes the expectation value using a sparse matrix representation of an operator for a density matrix.
Arguments
state::AbstractVectorS
: The quantum state vector.rho::sa.SparseMatrixCSC
: The density matrix.op::QuantumOps
: The quantum operator.op_str::String
: The string representation of the operator.qubit::Int
: The specific qubit index.matrix::sa.SparseMatrixCSC
: The sparse matrix representation of the operator.
Returns
- The expectation value as a
Float64
or a vector ofFloat64
for multiple qubits.
BlueTangle.correlation
— Functioncorrelation(m::Measurement, qubits::Vector) -> Any
Calculate correlation from a measurement on specified qubits.
m
: AMeasurement
object.qubits
: Vector of qubits for measuring correlations, e.g. eg: qubits=[1,3] measure <Z1Z3>
Returns the calculated correlation.
Alias:
correlation(state::AbstractVectorS, list_of_operators::String, qubits_applied::Vector) -> Float64
correlation(rho::sa.SparseMatrixCSC, list_of_operators::String, qubits_applied::Vector) -> Float64
Calculate the correlation for a given set of operators applied to specific qubits in either a quantum state vector or a density matrix. This function has two primary forms:
correlation(state::AbstractVectorS, list_of_operators::String, qubits_applied::Vector)
: Computes the correlation for a quantum state vector (state
) with a specified list of operators and qubits.correlation(rho::sa.SparseMatrixCSC, list_of_operators::String, qubits_applied::Vector)
: Computes the correlation for a density matrix (rho
) with a specified list of operators and qubits.
The corr_from_rho
function is an alias to get_corr
for density matrices.
Arguments
state::AbstractVectorS
: The quantum state vector.rho::sa.SparseMatrixCSC
: The density matrix.list_of_operators::String
: A string representing a list of operators, e.g., "Z,Z".qubits_applied::Vector
: A vector of qubit indices on which the operators are applied.
Returns
Float64
: The computed correlation value.
Examples
# For a state vector
state = sa.SparseVector([...]) # define your state vector
correlation = correlation(state, "Z,Z", [2, 4])
# For a density matrix
rho = sa.SparseMatrixCSC([...]) # define your density matrix
correlation = correlation(rho, "Z,Z", [2, 4])
Circuit functions
BlueTangle.compile
— Functioncompile(ops::Vector{<: QuantumOps}, options::Options=Options()) -> Circuit
Compile a set of quantum operations into a circuit.
ops
: A vector ofQuantumOps
objects.options
: Optional compilation options.
Returns a compiled Circuit
object.
BlueTangle.measure
— Functionmeasure(state::AbstractVectorS,number_of_experiment::Int)
this creates a measurement object from state vector.
BlueTangle.measure_ZNE
— Functionmeasure_ZNE(circuit::Circuit, number_of_experiment::Int) -> Measurement
Measure a quantum circuit multiple times.
circuit
: ACircuit
object.number_of_experiment
: Number of times to execute the circuit.
Returns a Measurement
object with the results.
BlueTangle.to_state
— Functionto_state(circuit::Circuit; init_state::AbstractVectorS=sa.sparse([])) -> sa.SparseVector
Convert a quantum circuit to a state vector.
circuit
: ACircuit
object.init_state
: (Optional) Initial state vector.
Returns a state vector representing the circuit.
to_state(MPS::it.MPS,M::Vector)
Convert MPS to state vector
BlueTangle.to_rho
— Functionto_rho(circuit::Circuit) -> sa.SparseMatrixCSC
Convert a quantum circuit to a density matrix (rho).
circuit
: ACircuit
object.
Returns a density matrix representing the circuit.
BlueTangle.to_MPS
— Functionto_MPS(vec::AbstractVector, sites::AbstractVector; kwargs...)
Convert state vector to MPS
Quantum Toolkit
BlueTangle.entanglement_entropy
— Functionentanglement_entropy(psi::AbstractVectorS) -> Float64
Calculates the entanglement entropy of a quantum state.
psi
: The quantum state vector (sa.SparseVector) for which to calculate the entropy.
Returns the calculated entanglement entropy.
entanglement_entropy(rho::sa.SparseMatrixCSC) -> Float64
Calculates the entanglement entropy of a density matrix.
rho
: The density matrix (sa.SparseMatrixCSC) for which to calculate the entropy.
Returns the calculated entanglement entropy.
BlueTangle.shadow
— Functionshadow(circuit::Circuit, number_of_experiment::Int) -> sa.SparseMatrixCSC
Construct a density matrix (rho) from classical shadow representation.
circuit
: ACircuit
object.number_of_experiment
: Number of experiments to run.
Returns a sparse density matrix representing the classical shadow.
BlueTangle.mag_moments
— Functionmag_moments(rho::sa.SparseMatrixCSC, op_str::String, moment_order::Int) -> Float64
Calculates the magnetization moments from a given density matrix.
rho
: The density matrix (sa.SparseMatrixCSC) of the quantum system.op_str
: String representation of the operator used for calculating the magnetization.moment_order
: The order of the magnetization moment to compute.
Returns the computed magnetization moment of the specified order.
calculates average magnetization moments from sample
Other
BlueTangle.get_N
— Functionget_N(state::AbstractVectorS) -> Int
get_N(rho::AbstractMatrixS) -> Int
qubit number from a set of operations
BlueTangle.fock_basis_create
— Functionfock_basis_create(N::Int)
BlueTangle.hilbert
— Functionhilbert(N::Int, mat::AbstractMatrix, qubit::Int, target_qubit::Int)
Constructs a sparse matrix representing the action of a quantum gate in a Hilbert space associated with a quantum system of N
qubits
The gate mat
is applied to the qubit
andtarget_qubit
. If qubit
is greater than target_qubit
, a controlled swap is performed before applying mat
.
Arguments
N::Int
: The number of qubits in the system.mat::AbstractMatrix
: The quantum gate to be applied.qubit::Int
: The qubit to which the gate is applied.target_qubit::Int
: The target qubit to which the gate is applied.
Returns
SparseMatrix
: The resulting sparse matrix representation of the gate operation.
hilbert(N::Int, mat::AbstractMatrix, qubit::Int)
Constructs a sparse matrix representing the action of a quantum gate in a Hilbert space associated with a quantum system of N
qubits.
Arguments
N::Int
: The number of qubits in the system.mat::AbstractMatrix
: The quantum gate to be applied.qubit::Int
: The qubit to which the gate is applied.
Returns
SparseMatrix
: The resulting sparse matrix representation of the gate operation.
BlueTangle.expand_multi_op
— Functionexpand_multi_op(list_of_operators::String, qubits_applied::Vector, N::Int) -> Matrix
Expand multiple quantum operators over a specified set of qubits.
list_of_operators
: A string representing a list of operators.qubits_applied
: A vector of qubits on which operators are applied.N
: Total number of qubits.
Returns a matrix representing the expanded operators.
BlueTangle.string_to_matrix
— Functionstring_to_matrix(list_of_operators::String) -> Matrix
Convert a string of comma-separated operators into a matrix representation.
list_of_operators
: A string representing a list of operators, e.g.: "Z,Z,sa.IP(.2)"
Returns a matrix representation of the operators.
BlueTangle.zero_state
— Functionzero_state(N::Int) -> sa.SparseVector
Returns a sparse vector representing the |000...> quantum state.
create all zero state
BlueTangle.one_state
— Functioncreate all one state
BlueTangle.product_state
— Functionproduct_state(list_of_qubits::Vector) -> sa.SparseVector
Creates a quantum state vector from a list of qubits.
list_of_qubits
: A vector representing the state of each qubit.
Returns a sparse vector representing the quantum state of the system.
create given product state
BlueTangle.neel_state01
— Functioncreate neel state 010101
BlueTangle.neel_state10
— Functioncreate neel state 101010
BlueTangle.apply
— Functionapply(state::AbstractVectorS, op::QuantumOps)
Apply a quantum gate operation to a state vector in place.
state
: A sparse quantum state vector to be modified.op
: AQuantumOps
object representing the gate operation.
Modifies the state vector directly.
apply(rho::sa.SparseMatrixCSC, op::QuantumOps)
Apply a quantum gate operation to a state vector in place.
rho
: A sparse quantum density matrix to be modified.op
: AQuantumOps
object representing the gate operation.
Modifies the state vector directly.
BlueTangle.linear_fit
— Functionlinear_fit(xdata::Vector, ydata::Vector)
Compute the coefficients a and b for the linear fit of the given data.
Arguments
xdata
: Array of x values.ydata
: Array of corresponding y values.
Returns
a
: Intercept of the linear fit.b
: Slope of the linear fit.se_a
: Standard error of the intercept.se_b
: Standard error of the slope.
Description
This function fits a simple linear model (y = a + b*x) to the provided data points. It checks if the lengths of xdata and ydata are the same and then calculates the coefficients for the linear fit. Additionally, it computes the standard errors for both the slope and the intercept.
The function uses the least squares method for the linear regression. The standard errors are calculated based on the residual sum of squares and the total sum of squares for the x values.
BlueTangle.quadratic_fit
— Functionquadratic_fit(xdata::Vector{Float64}, ydata::Vector{Float64}) -> Vector{Float64}
Fit a quadratic function to the given data.
Arguments
xdata::Vector{Float64}
: A vector of x-coordinates.ydata::Vector{Float64}
: A vector of y-coordinates corresponding toxdata
.
Returns
Vector{Float64}
: The coefficients[a, b, c]
of the fitted quadratic functiony = ax^2 + bx + c
.
Description
This function performs a least squares quadratic fit to the input data.
Object types
BlueTangle.QuantumOps
— TypeQuantumOps
Abstract type representing quantum operations.
BlueTangle.Op
— TypeOp(q::Int, name::String, mat::AbstractMatrix, qubit::Int, target_qubit::Int, noise::QuantumChannel) <: QuantumOps
Represents a quantum operation.
q
: Number of qubits involved in the operation.name
: Name of the operation.mat
: Matrix representation of the quantum operation.qubit
: index of the target qubit.target_qubit
: index of the target qubit for two-qubit operations.noise
: Noise model associated with the operation.
Constructs an Op object representing a quantum operation with optional noise.
BlueTangle.ifOp
— TypeifOp(q::Int, name::String, mat::AbstractMatrix, qubit::Int, if01::Tuple{Matrix,Matrix}, noise::QuantumChannel) <: QuantumOps
Represents a conditional quantum operation used for mid-circuit born measurements. It is specifically designed for mid-circuit measurements in the X, Y, Z, or a random basis (R). Depending on the measurement outcomes (0 or 1), different gates specified in if01
can be applied.
Fields
q
: Integer representing the number of qubits involved in the operation.name
: String indicating the name of the operation. Valid names are "MX", "MY", "MZ", or "MR", corresponding to operations in the X, Y, Z basis, or a random basis (R), respectively.mat
: Matrix representing the quantum operation.qubit
: Integer specifying the index of the target qubit.if01
: Tuple of two matrices. The first matrix is applied if the measured state of the qubit is 0, and the second matrix is applied if the measured state is 1.noise
: Represents a 'born measurement quantum channel', indicating the noise model associated with the operation.
Usage
ifOp
is constructed to represent quantum operations that are conditional on the measurement outcome of a qubit. This operator is particularly useful in quantum circuits for implementing dynamic responses based on mid-circuit measurement results.
Examples
```julia
Example of creating an ifOp for a conditional operation based on the X-basis measurement (MX)
conditionalop = ifOp("MX", qubit, (operationif0, operationif_1))
BlueTangle.Measurement
— TypeMeasurement(bitstr::Union{Vector, UnitRange}, sample::Vector, expect::Vector, mag_moments::Vector, measurement_basis::String, number_of_experiment::Int, circuit_name::String, number_of_qubits::Int, density_matrix::sa.SparseMatrixCSC)
Represents the result of quantum measurements.
bitstr
: Basis of measurement represented as integers or a range.sample
: Vector of probabilities.expect
: Expectation values of the measurement.mag_moments
: Magnetic moments obtained from the measurement.measurement_basis
: Measurement basis used.number_of_experiment
: Number of experiments performed.circuit_name
: Name of the circuit used.number_of_qubits
: Number of qubits involved.density_matrix
: Density matrix obtained from the measurement.
Constructs a Measurement object to store the results of quantum measurements.
BlueTangle.QuantumChannel
— TypeQuantumChannel(q::Int, model::String, p::Float64)
Represents a quantum noise channel.
q
: Number of qubits affected by the noise channel.model
: The noise model type as a string.p
: Probability parameter for the noise model.kraus
: A vector of matrices representing Kraus operators for the channel.
Constructs a QuantumChannel object for specified qubits, noise model, and probability.
BlueTangle.NoiseModel
— TypeCreate a NoiseModel
BlueTangle.AnsatzOptions
— TypeAnsatzOptions(; N::Int, ops::Union{Vector{String},Vector{<:QuantumOps}}, noise=false, init::Union{sa.SparseVector,Circuit}=sa.sparse([]), model::String="lbfgs", number_of_iterations::Int=1000, learning_rate::Float64=0.01, pars_initial::Vector=[], deep_circuit::Bool=false, history::Bool=true)
Constructs an AnsatzOptions
object that contains the configuration for a variational quantum circuit ansatz.
Arguments
N::Int
: The number of qubits in the ansatz.ops::Union{Vector{String},Vector{<:QuantumOps}}
: The quantum operations defining the ansatz, either as a vector of strings or a vector ofQuantumOps
.noise=false
: Specifies whether to include noise in the ansatz. Can be either aNoiseModel
or a boolean value.init::Union{sa.SparseVector,Circuit}=sa.sparse([])
: The initial state of the ansatz, either as a sparse vector or aCircuit
. Defaults to the zero state.
Optimization Models
Gradient-Free Models (uses PRIMA.jl)
"cobyla"
: Affine model, supports bounds, linear, and non-linear constraints."newuoa"
: Quadratic model, no constraints."uobyqa"
: Quadratic model, no constraints."bobyqa"
: Quadratic model, supports bounds."lincoa"
: Quadratic model, supports bounds and linear constraints."prima"
: Automatically picks the best PRIMA model.
Gradient-Based Models (uses ForwardDiff.jl, Optimisers.jl, OptimKit.jl]
"lbfgs"
: Uses LBFGS optimizer (gradient-based)."adam"
: Uses Adam optimizer (gradient-based)."descent"
: Uses gradient descent optimizer."radam"
: Uses RAdam optimizer."momentum"
: Uses Momentum optimizer."nesterov"
: Uses Nesterov optimizer.number_of_iterations::Int=1000
: The number of iterations for the optimization.learning_rate::Float64=0.01
: The learning rate for the optimization.pars_initial::Vector=[]
: The initial parameters for the ansatz. If not provided, random parameters are generated.deep_circuit::Bool=false
: Specifies whether to use a deep circuit ansatz.history::Bool=true
: Specifies whether to record the optimization history.
Returns
- An
AnsatzOptions
object containing the configuration for the variational quantum circuit ansatz.
BlueTangle.Circuit
— TypeCircuit
Represents a quantum circuit.
stats
: NamedTuple containing statistics about the circuit.options
: Options object with circuit configurations.ops
: Vector of QuantumOps representing the operations in the circuit.
Constructs a Circuit object representing a quantum circuit.
BlueTangle.Options
— TypeOptions(circuit_name::String, measurement_basis::String, final_measurement_noise::QuantumChannel, Noise1::QuantumChannel, Noise2::QuantumChannel, twirl::Bool, noisy_swap::Bool, density_matrix::Bool)
Represents configuration options for a quantum circuit.
circuit_name
: Name of the circuit.measurement_basis
: Measurement basis used in the circuit.final_measurement_noise
: Noise model for final measurement error.noise
: Noise model for single-qubit and two-qubit operations.twirl
: Boolean flag for twirling operations.noisy_swap
: Boolean flag for swap errors.density_matrix
: Boolean flag to indicate if density matrix should be calculated.
Constructs an Options object with specified settings for a quantum circuit.
BlueTangle.Layout
— TypeCreate Layout
Gates
BlueTangle.gate
— Constantgate
A constant that holds common quantum gates and projectors.
This includes:
- Identity (
I
) - Pauli gates (
X
,Y
,Z
) - Hadamard (
H
) - Phase gates (
S
,T
) - Special gates like
sqrt(X)
(equal to RX(pi/2)exp(1impi/4)) - Projectors (
P0
for |0><0|,P1
for |1><1|) - Controlled gates such as
CX
(CNOT),CNOT
(an alias for CX), andCZ
- The
SWAP
,ISWAP
,FSWAP
gate - The
ECR
,SYC
gate
Each single-qubit gate is represented as a 2x2 matrix, while multi-qubit gates like CNOT
, ECR
, SYC
, CZ
, and SWAP
, ISWAP
, FSWAP
are represented as 4x4 matrices.
BlueTangle.random_ops
— Functionrandom_ops_2(N::Int, len::Int; measure_prob::Float64=0.0, measure_basis::Vector{String}=["MX", "MY", "MZ"]) -> Vector{QuantumOps}
Create a sequence of random quantum gate operations, with optional mid-circuit measurements.
Arguments
N::Int
: The number of qubits in the system.len::Int
: The length of the sequence of operations to generate.
Keyword Arguments
measure_prob::Float64
: The probability of adding a measurement operation after each gate.measure_basis::Vector{String}
: The basis in which measurements are performed.
Returns
Vector{QuantumOps}
: A vector of randomly chosen quantum operations (QuantumOps
), each representing a gate or a measurement operation.
Description
This function creates a vector of quantum operations, where each operation is either a randomly chosen gate from the set {"X", "Y", "Z", "H", "S", "CX","CZ","CP","GIVENS","FSIM","SWAP","ISWAP","FSWAP","SYC","ECR"} or a measurement operation, based on measure_prob
.
Example
ops = random_ops(5, 10; measure_prob=0.2, measure_basis=["MX","MZ"])
This example generates a sequence of 10 random gates and measurements (with a 20% chance of a measurement after each gate) for a 5-qubit system.
BlueTangle.random_clifford
— Functionrandom_clifford(N::Int, len::Int; measure_prob::Float64=0.0, measure_basis::Vector{String}=["MX","MY","MZ"]) -> Vector{QuantumOps}
Generate a random sequence of Clifford gates, with optional mid-circuit measurements.
Arguments
N::Int
: The number of qubits in the system.len::Int
: The length of the sequence of operations to generate.
Keyword Arguments
measure_prob::Float64
: The probability of adding a measurement operation after each gate.measure_basis::Vector{String}
: The basis in which measurements are performed.
Returns
Vector{QuantumOps}
: A vector of randomly chosen quantum operations (QuantumOps
), each representing a Clifford gate or a measurement operation.
Description
This function creates a vector of quantum operations, where each operation is either a randomly chosen Clifford gate from the set {"X", "Y", "Z", "H", "S", "CNOT", "SWAP", "CZ"} or a measurement operation, based on measure_prob
. For two-qubit gates ("CNOT", "SWAP", "CZ"), adjacent qubits (qubit r
and qubit r+1
) are selected. For single-qubit gates and measurements, a random qubit r
is chosen.
Example
ops = random_clifford(5, 10; measure_prob=0.2, measure_basis=["MX","MZ"])
This example generates a sequence of 10 random Clifford gates and measurements (with a 20% chance of a measurement after each gate) for a 5-qubit system.
BlueTangle.get_stats
— Functionget statistics about the operations
BlueTangle.get_layers
— Functionget_layers for a given set of ops
Noise
BlueTangle.is_valid_quantum_channel
— Functionis_valid_quantum_channel(kraus::Vector{Matrix}) -> Bool
Determines the validity of a set of Kraus operators.
kraus
: A vector of matrices, each representing a Kraus operator.
This function checks if the provided Kraus operators form a valid quantum channel. It does so by verifying if the sum of the products of each Kraus operator and its adjoint (approximately) equals the identity matrix. Returns true
if the set is valid, false
otherwise.
BlueTangle.Noise1
— FunctionNoise1(model::String, p::Float64)
Constructors for QuantumChannel objects for 1 qubit.
model
: The noise model type as a string.p
: Probability parameter for the noise model.
Returns a QuantumChannel object.
BlueTangle.Noise2
— FunctionNoise2(model::String, p::Float64)
Constructors for QuantumChannel objects for 2 qubits.
model
: The noise model type as a string.p
: Probability parameter for the noise model.
Returns a QuantumChannel object.
BlueTangle.apply_noise
— Functionapply noise on qubit or target_qubit of a given state and noise model
apply noise on qubit or target_qubit of a given density matrix and noise model
BlueTangle.apply_twirl
— Functionapply_twirl(op::QuantumOps) -> Vector{QuantumOps}
Applies twirling with specified noise to a single quantum operation.
op
: A QuantumOps object representing the quantum operation.
Returns a vector of QuantumOps representing the twirled operation with noise.
BlueTangle.cnot_amplifier!
— Functioncnot_amplifier!(ops::Vector{<:QuantumOps}, CNOT_pair=0)
Amplify the presence of CNOT (or CX) operations in a vector of quantum operations.
This function adds extra pair of CNOT operation in the ops
vector a specific number of times in place. This is useful for amplifying the effect of noise via CNOT operations in a sequence of quantum operations.
Arguments
ops::Vector{<:QuantumOps}
: A vector of quantum operations, whereT
is a subtype ofQuantumOps
.CNOT_pair::Int
(optional): The number of additional pairs of CNOT operations to insert for each original CNOT operation inops
. The default value is 0, which means no additional operations are inserted.
Examples
ops = [Op("H",1), Op("CNOT",1,2), Op("X",2)]
cnot_amplifier!(ops, 1)
# `ops` will be modified to: [Op("H",1), Op("CNOT",1,2), Op("CNOT",1,2), Op("CNOT",1,2), Op("X",2)]
BlueTangle.op_amplifier!
— Functionop_amplifier!(ops::Vector{<:QuantumOps},op_pair=0)
same as cnot_amplifier!
but for all operations
BlueTangle.error_mitigate_data
— Functionerrormitigatedata(xdata::Vector, ydata::Vector)
Perform error mitigation on a dataset by fitting a linear model and extracting the estimate and standard error.
This function takes two vectors xdata
and ydata
which represent the independent and dependent variables of a dataset, respectively.
Arguments
xdata::Vector
: The independent variable data points.ydata::Vector
: The dependent variable data points, corresponding to each xdata point.
Returns
est
: The estimated intercept from the linear fit.se
: The standard error of the estimated intercept.fit_plot
: A tuple containing the x-values from 0 to the last element ofxdata
and the corresponding fitted y-values from the model.
Hamiltonian
BlueTangle.hamiltonian
— Functionhamiltonian(N::Int, string_of_ops::Vector, boundary::String="open")
Constructs a Hamiltonian matrix for a quantum system with N
qubits. The Hamiltonian is built based on the operators and their corresponding couplings specified in string_of_ops
.
The string_of_ops
should be an alternating array of coupling constants and operator strings. For example, [.1, "Z,Z", .5, "X"]
implies a system with alternating couplings .1
and .5
, and operators "Z,Z" and "X".
The boundary
parameter specifies the boundary conditions of the system. It can be either "open" or "periodic". In the case of "open" boundary conditions, interactions are only included for sites within the system size. For "periodic" boundary conditions, the system is treated as if it were wrapped around itself, allowing interactions that cross the end and start of the chain.
Arguments
N::Int
: The number of qubits in the system.string_of_ops::Vector
: An alternating vector of coupling constants and operator strings.boundary::String
: The boundary condition of the system, either "open" or "periodic" (default is "open").
Returns
SparseMatrix
: The Hamiltonian matrix representing the specified quantum system.
Example
N = 6
string_of_ops = [-J, "Z,Z", -h, "X"]
H = hamiltonian(N, string_of_ops, "open")
This function iterates over the operator strings, applies each operator to the appropriate qubits based on the boundary conditions, and scales them by their corresponding coupling constants to construct the Hamiltonian.
hamiltonian(rows_cols::Union{Tuple{Int64, Int64},Vector{Int64}}, string_of_ops::Vector, boundary::String="open")
Constructs a 2D Hamiltonian matrix for a quantum system with rows * cols
qubits, potentially with double counting.
Arguments
rows_cols
: A tuple or vector specifying the dimensions of the 2D lattice (rows, columns).string_of_ops::Vector
: An alternating vector of coupling constants and operator strings.boundary::String
: The boundary condition of the system, either "open" or "periodic" (default is "open").
Returns
SparseMatrixCSC
: The Hamiltonian matrix representing the specified 2D quantum system.
BlueTangle.hamiltonian_exp
— Functionhamiltonian_exp(N::Int, total_time::Float64, string_of_ops::Vector; dt=0.1) -> Vector{QuantumOps}
Expands a given Hamiltonian expressed as a string of operations into a sequence of quantum gates using Trotterization.
N
: Number of qubits involved in the simulation.total_time
: The total simulation time over which the Hamiltonian is to be applied.string_of_ops
: A vector where odd indices contain the coupling strengths and even indices contain comma-separated strings representing the operators (X, Y, Z) applied to consecutive qubits.dt
(optional): The time step for Trotterization, with a default value of 0.1.
The function hamiltonian_exp
parses the string_of_ops
to construct a sequence of operations based on the specified operators and their coupling strengths. For each term in the Hamiltonian:
This method returns a vector of QuantumOps
, each representing a quantum operation to be applied sequentially to simulate the Hamiltonian over the specified time total_time
.
Example
Define a Hamiltonian for a 3-qubit system with mixed interactions over 2 seconds
N = 3 totaltime = 2.0 stringofops = [1.0, "X,Y", 0.5, "Y,Z"] ops = hamiltonianexp(N, totaltime, stringof_ops)
ops will contain a sequence of quantum gates to apply.
hamiltonian_exp(rows_cols::Union{Tuple{Int64, Int64}, Vector{Int64}}, total_time::Float64, string_of_ops::Vector;
dt=0.01, full=true, enumeration_style::Symbol=:rowwise) -> Vector{QuantumOps}
Generate a sequence of quantum operations to simulate a 2D lattice Hamiltonian using Trotterization.
Arguments
rows_cols
: A tuple or vector specifying the dimensions of the 2D lattice (rows, columns).total_time
: The total simulation time.string_of_ops
: A vector alternating between coupling strengths and operator strings. Each operator string is a comma-separated list of Pauli operators (X, Y, Z) representing the terms in the Hamiltonian.
Keyword Arguments
dt
: Time step for Trotterization (default: 0.01).full
: If true, repeats the operation sequence for each Trotter step (default: true).enumeration_style
: Specifies how qubits are indexed in the 2D lattice. Options are:rowwise
(default) or:colwise
.
Returns
A vector of QuantumOps
representing the sequence of quantum operations for the simulation.
Details
This function implements a Trotterized evolution of a 2D lattice Hamiltonian. It supports arbitrary Pauli string operators and handles both horizontal and vertical nearest-neighbor interactions in the lattice.
The Hamiltonian terms are applied to all relevant qubit pairs in the lattice. For operators involving two qubits (e.g., "X,Y"), both horizontal and vertical applications are considered. Single-qubit operators are applied to each qubit individually.
The function uses the _apply_term
helper function to generate the specific quantum operations for each term in the Hamiltonian.
Example
rows, cols = 3, 4
total_time = 2.0
string_of_ops = [1.0, "X,Y", 0.5, "Z"]
ops = hamiltonian_exp((rows, cols), total_time, string_of_ops)
BlueTangle.VQE
— FunctionVQE(opt::AnsatzOptions)
Performs the Variational Quantum Eigensolver (VQE) algorithm to find the optimal state of a given loss function using the variational quantum circuit defined by the AnsatzOptions
.
Arguments
opt::AnsatzOptions
: TheAnsatzOptions
object containing the configuration for the variational quantum circuit.
Returns
- A tuple containing:
- The optimization history (energy values) if
history=true
, otherwise the final energy value. - The optimized parameters for the variational circuit.
- The final state obtained from the optimized parameters.
- The optimization history (energy values) if
BlueTangle.variational_apply
— Functionvariational_apply(pars::AbstractVectorS, opt::AnsatzOptions)
Applies the variational quantum circuit defined by the AnsatzOptions
to the initial state using the given parameters.
Arguments
pars::Vector
: The parameters for the variational circuit.opt::AnsatzOptions
: TheAnsatzOptions
object containing the configuration for the variational circuit.
Returns
- The final state after applying the variational circuit.
Index
BlueTangle.gate
BlueTangle.AnsatzOptions
BlueTangle.Circuit
BlueTangle.Layout
BlueTangle.Measurement
BlueTangle.NoiseModel
BlueTangle.Op
BlueTangle.Options
BlueTangle.QuantumChannel
BlueTangle.QuantumOps
BlueTangle.ifOp
BlueTangle.Noise1
BlueTangle.Noise2
BlueTangle.VQE
BlueTangle.apply
BlueTangle.apply_noise
BlueTangle.apply_twirl
BlueTangle.cnot_amplifier!
BlueTangle.compile
BlueTangle.correlation
BlueTangle.entanglement_entropy
BlueTangle.error_mitigate_data
BlueTangle.expand_multi_op
BlueTangle.expect
BlueTangle.fock_basis_create
BlueTangle.get_N
BlueTangle.get_layers
BlueTangle.get_probs_from_sample
BlueTangle.get_stats
BlueTangle.hamiltonian
BlueTangle.hamiltonian_exp
BlueTangle.hilbert
BlueTangle.is_valid_quantum_channel
BlueTangle.linear_fit
BlueTangle.mag_moments
BlueTangle.measure
BlueTangle.measure_ZNE
BlueTangle.neel_state01
BlueTangle.neel_state10
BlueTangle.one_state
BlueTangle.op_amplifier!
BlueTangle.product_state
BlueTangle.quadratic_fit
BlueTangle.random_clifford
BlueTangle.random_ops
BlueTangle.sample
BlueTangle.shadow
BlueTangle.string_to_matrix
BlueTangle.to_MPS
BlueTangle.to_rho
BlueTangle.to_state
BlueTangle.variational_apply
BlueTangle.zero_state