QuantumUtilities

Documentation for QuantumUtilities.

QuantumUtilities.BosonicEnvironmentMethod
BosonicEnvironment(J, kT)

Construct a bosonic environment with spectral density J and temperature kT.

Arguments

  • J: The spectral density. Must be a subtype of AbstractSD.
  • kT: The temperature (times the Boltzmann constant).

Returns

  • The bosonic environment.
source
QuantumUtilities.EnvironmentMethod
Environment{E}(J, kT, μ)

Construct an environment of type E, with spectral density J, temperature kT, and chemical potential μ.

Arguments

  • J: The spectral density. Must be a subtype of AbstractSD.
  • kT: The temperature (times the Boltzmann constant).
  • μ: The chemical potential.

Returns

  • The environment.
source
QuantumUtilities.EnvironmentMethod
Environment(E, J, kT, μ)

Construct an environment of type E, with spectral density J, temperature kT, and chemical potential μ.

Arguments

  • E: Type of the environment. Must be as subtype of EnvironmentType.
  • J: The spectral density. Must be a subtype of AbstractSD.
  • kT: The temperature (times the Boltzmann constant).
  • μ: The chemical potential.

Returns

  • The environment.
source
QuantumUtilities.EnvironmentMethod
Environment{E}(J, kT)

Construct an environment of type E, with spectral density J and temperature kT.

Arguments

  • J: The spectral density. Must be a subtype of AbstractSD.
  • kT: The temperature (times the Boltzmann constant).

Returns

  • The environment.
source
QuantumUtilities.EnvironmentMethod
Environment(E, J, kT)

Construct an environment of type E, with spectral density J and temperature kT.

Arguments

  • E: Type of the environment. Must be as subtype of EnvironmentType.
  • J: The spectral density. Must be a subtype of AbstractSD.
  • kT: The temperature (times the Boltzmann constant).

Returns

  • The environment.
source
QuantumUtilities.FermionicEnvironmentMethod
FermionicEnvironment(J, kT, μ)

Construct a bosonic environment with spectral density J, temperature kT, and chemical potential μ.

Arguments

  • J: The spectral density. Must be a subtype of AbstractSD.
  • kT: The temperature (times the Boltzmann constant).
  • μ: The chemical potential.

Returns

  • The fermionic environment.
source
QuantumUtilities.SpinLengthType
struct SpinLength{N,D}

A type representing a spin length N/D.

Note that N and D must be such that the ratio N/D is either an integer or a half-integer.

source
QuantumUtilities.SpinLengthMethod
SpinLength(N, D)
SpinLength(S0::Int)
SpinLength(S0::Rational)
SpinLength(S0)
SpinInteger(S0)
SpinHalfInteger(S0)

Convenient constructors for spin length types.

source
QuantumUtilities.addMethod
add(S1::SpinLength, S2::SpinLength)

Return a tuple of the different possible total angular momentum values resulting from the addition of two angular momenta S1 and S2.

Arguments

  • S1::SpinLength: The first angular momentum.
  • S2::SpinLength: The second angular momentum.

Returns

  • A tuple containing all possible values of total angular momentum.
source
QuantumUtilities.annihilation_operatorFunction
annihilation_operator(ncutoff::Int, noffset::Int=0)

Create the annihilation (lowering) operator matrix in the Fock basis for a quantum harmonic oscillator.

Arguments

  • ncutoff::Int: The cutoff value for the number of excitations to be included.
  • noffset::Int=0: The offset for the number of excitations, defaults to 0.

Returns

  • The matrix representation of the annihilation operator.
source
QuantumUtilities.anticommutator_superopMethod
anticommutator_superop(A)

Computes the Liouville space superoperator representation of the anti-commutator with a matrix A.

Arguments

  • A: The input matrix.

Returns

The anti-commutator superoperator matrix obtained by adding the right superoperator of A to the left superoperator of A. ```

source
QuantumUtilities.cauchy_quadgkMethod
cauchy_quadgk(g, a, b; kws...)

Computes the Cauchy principal value of the integral of a function g over the interval [a, b] using the quadgk quadrature method.

Arguments

  • g: The function to integrate.
  • a: The lower bound of the interval.
  • b: The upper bound of the interval.
  • kws...: Additional keyword arguments accepted by quadgk.

Returns

A tuple (I, E) containing the approximated integral I and an estimated upper bound on the absolute error E.

Throws

  • ArgumentError: If the interval [a, b] does not include zero.

Examples

julia> g(x) = 1 / (x + 1)

julia> cauchy_quadgk(g, -1/2, 1/2)
(-1.09861228866811, 1.9939216944209193e-11)
source
QuantumUtilities.coherent_stateFunction
coherent_state(α, ncutoff::Int, noffset::Int=0)

Create a coherent state α in the Fock basis of a quantum harmonic oscillator.

The state is created by applying the displacement operator to the vacuum state:

\[\left|\alpha\right\rangle = D(\alpha)\left|0\right\rangle.\]

Arguments

  • α: The complex amplitude of the coherent state.
  • ncutoff::Int: The cutoff value for the number of excitations to be included.
  • noffset::Int=0: The offset for the number of excitations, defaults to 0.

Returns

  • The coherent state.
source
QuantumUtilities.coherent_state_analyticFunction
coherent_state_analytic(α, ncutoff::Int, noffset::Int=0)

Create a coherent state α in the Fock basis of a quantum harmonic oscillator.

The state is created by using the analytical expression for the state coefficients:

\[\left|\alpha\right\rangle = \sum_{n=\mathrm{noffset}}^{\mathrm{ncutoff}} e^{-\frac{|\alpha|^2}{2}}\frac{\alpha^n}{\sqrt{n!}}\left|n\right\rangle.\]

Due to the truncation of the series, the state is not guaranteed to be normalized.

Arguments

  • α: The complex amplitude of the coherent state.
  • ncutoff::Int: The cutoff value for the number of excitations to be included.
  • noffset::Int=0: The offset for the number of excitations, defaults to 0.

Returns

  • The coherent state.
source
QuantumUtilities.coherent_state_coefficientMethod
coherent_state_coefficient(α, n::Int)

Coefficient in the Fock basis for the coherent state α of a quantum harmonic oscillator:

\[\left\langle n\middle|\alpha\right\rangle = e^{-\frac{|\alpha|^2}{2}}\frac{\alpha^n}{\sqrt{n!}}.\]

Arguments

  • α: The complex amplitude of the coherent state.
  • n::Int: The Fock basis occupation number of the coefficient.

Returns

  • The coherent state coefficient in the Fock basis.
source
QuantumUtilities.commutator_superopMethod
commutator_superop(A)

Computes the Liouville space superoperator representation of the commutator with a matrix A.

Arguments

  • A: The input matrix.

Returns

The commutator superoperator matrix obtained by subtracting the right superoperator of A from the left superoperator of A. ```

source
QuantumUtilities.creation_operatorFunction
creation_operator(ncutoff::Int, noffset::Int=0)

Create the creation (raising) operator matrix in the Focks basis for a quantum harmonic oscillator.

Arguments

  • ncutoff::Int: The cutoff value for the number of excitations to be included.
  • noffset::Int=0: The offset for the number of excitations, defaults to 0.

Returns

  • The matrix representation of the creation operator.
source
QuantumUtilities.displacement_operatorFunction
displacement_operator(α, ncutoff::Int, noffset::Int=0)

Create the displacement operator for a quantum harmonic oscillator.

The operator is constructed using the definition as unitary operator based on the creation and annihilation operators:

\[D(\alpha) = e^{\alpha a^\dagger - \alpha^* a}.\]

Arguments

  • α: The complex displacement parameter.
  • ncutoff::Int: The cutoff value for the number of excitations to be included.
  • noffset::Int=0: The offset for the number of excitations, defaults to 0.

Returns

  • The matrix representation of the displacement operator.
source
QuantumUtilities.displacement_operator_analyticFunction
displacement_operator_analytic(α, ncutoff::Int, noffset::Int=0)

Create the displacement operator for a quantum harmonic oscillator.

The operator is constructed using the analytic expressions of its matrix elements in the Fock basis:

\[\left\langle m \middle| D(\alpha) \middle| n \right\rangle = e^{-\frac{|\alpha|^2}{2}} \sqrt{\frac{n!}{m!}} \alpha^{m-n} L_n^{(m-n)}(|\alpha|^2),\]

where $L_n^{(m-n)}$ are the generalised Laguerre polynomial.

Due to the truncation of the matrix elements, the operator is not guaranteed to be unitary.

Arguments

  • α: The complex displacement parameter.
  • ncutoff::Int: The cutoff value for the number of excitations to be included.
  • noffset::Int=0: The offset for the number of excitations, defaults to 0.

Returns

  • The matrix representation of the displacement operator.
source
QuantumUtilities.displacement_operator_coefficientMethod
displacement_operator_coefficient(α, m::Int, n::Int)

Return the displacement operator matrix elements in the Fock basis for a quantum harmonic oscillator:

\[\left\langle m \middle| D(\alpha) \middle| n \right\rangle = e^{-\frac{|\alpha|^2}{2}} \sqrt{\frac{n!}{m!}} \alpha^{m-n} L_n^{(m-n)}(|\alpha|^2),\]

where math L_n^{(m-n)} are the generalised Laguerre polynomial.

Arguments

  • α: The complex displacement parameter.
  • m::Int: The occupation number of the row.
  • n::Int: The occupation number of the column.

Returns

  • The matrix elements of the displacement operator in the Fock basis.
source
QuantumUtilities.dissipator_superopMethod
dissipator_superop(L, M)
dissipator_superop(L)

Computes the superoperator corresponding to the Lindblad dissipator

\[D(\cdot) = L \cdot M^\dagger - \frac{1}{2} \{M^\dagger L, \cdot\}\]

If only L is provided, it is assumed that M = L.

Arguments

  • L: The jump operator.

Returns

The superoperator matrix representing the Lindblad dissipator. ```

source
QuantumUtilities.hamiltonian_evolution_superopMethod
hamiltonian_evolution_superop(H)

Computes the superoperator corresponding to the Hamiltonian contribution to the time evolution given by the Hamiltonian matrix H

\[L(\cdot) = -i[H, \cdot]\]

Arguments

  • H: The Hamiltonian matrix representing the dynamics of the system.

Returns

The superoperator matrix representing the Hamiltonian evolution. ```

source
QuantumUtilities.left_right_superopMethod
left_right_superop(A, B)

Computes the Liouville space superoperator representation of the operation A ⋅ B.

Arguments

  • A: The left matrix.
  • B: The right matrix.

Returns

The superoperator obtained by performing a Kronecker product between the transpose of A and B. ```

source
QuantumUtilities.left_superopMethod
left_superop(A)
left_superop(A, d::Int)

Computes the Liouville space left superoperator representation of a matrix A.

Arguments

  • A: The input matrix.
  • d: (Optional) The size of the identity matrix. If not provided, it is set to the number of columns in A.

Returns

The left superoperator matrix obtained by performing a Kronecker product between the identity matrix and A.

Examples

julia> A = [1 2; 3 4]
2×2 Matrix{Int64}:
 1  2
 3  4

julia> left_superop(A)
4×4 Matrix{Int64}:
 1  2  0  0
 3  4  0  0
 0  0  1  2
 0  0  3  4
source
QuantumUtilities.lindbladian_superopMethod
lindbladian_superop(H, Ls, γs::AbstractVector)
lindbladian_superop(H, Ls, γs::AbstractMatrix)
lindbladian_superop(H)

Computes the Lindblad superoperator for the given Hamiltonian H, jump operators Ls, and decay rates γs.

\[L = -i[H, \cdot] + \sum_i γ_{n,m} (L_n \cdot L_m^\dagger - \frac{1}{2} \{L_m^\dagger L_n, \cdot\})\]

If only H is provided, the Lindblad superoperator is equivalent to the free Hamiltonian evolution superoperator.

Arguments

  • H: The Hamiltonian.
  • Ls: The jump operators.
  • γs: The decay rates.

Returns

The superoperator matrix representing the Lindblad dissipator. ```

source
QuantumUtilities.momentum_operatorFunction
momentum_operator(ncutoff::Int, noffset::Int=0)

Create the momentum operator matrix in the Fock basis for a quantum harmonic oscillator.

Arguments

  • ncutoff::Int: The cutoff value for the number of excitations to be included.
  • noffset::Int=0: The offset for the number of excitations, defaults to 0.

Returns

  • The matrix representation of the momentum operator.
source
QuantumUtilities.number_operatorFunction
number_operator(ncutoff::Int, noffset::Int=0)

Create the number operator matrix in the Fock basis for a quantum harmonic oscillator.

Arguments

  • ncutoff::Int: The cutoff value for the number of excitations to be included.
  • noffset::Int=0: The offset for the number of excitations, defaults to 0.

Returns

  • The matrix representation of the number operator.
source
QuantumUtilities.occupationMethod
occupation(env::BosonicEnvironment, ω)

Return the occupation number of the bosonic environment env at frequency ω, that is

\[n_\mathrm{B}(\omega) = \frac{1}{e^{\frac{\omega}{\mathrm{kT}}} - 1},\]

where \mathrm{kT} is the temperature of env.

Arguments

  • env: The bosonic environment.
  • ω: The frequency.

Returns

  • The bosonic occupation at the given frequency.
source
QuantumUtilities.occupationMethod
occupation(env::FermionicEnvironment, ω)

Return the occupation number of the fermionic environment env at frequency ω, that is

\[n_\mathrm{F}(\omega) = \frac{1}{e^{\frac{\omega-\mu}{\mathrm{kT}}} + 1},\]

where \mathrm{kT} is the temperature and \mu the chemical potential of env.

Arguments

  • env: The fermionic environment.
  • ω: The frequency.

Returns

  • The fermionic occupation at the given frequency.
source
QuantumUtilities.operator_to_vectorMethod
operator_to_vector(A)

Converts a matrix or array-like object A into a one-dimensional vector.

Arguments

  • A: The input matrix or array-like object.

Returns

A one-dimensional vector representing the elements of A.

Examples

julia> A = [1 2; 3 4]
2×2 Matrix{Int64}:
 1  2
 3  4

julia> operator_to_vector(A)
4-element Vector{Int64}:
 1
 3
 2
 4
source
QuantumUtilities.partial_traceMethod
partial_trace(ρ::AbstractMatrix, trace_indices::Union{Int,Tuple{Vararg{Int}}}, dims::Tuple{Vararg{Int}})

Computes the partial trace of a density matrix ρ by tracing out the specified subsystems.

Arguments

  • ρ::AbstractMatrix: The input density matrix.
  • trace_indices::Union{Int,Tuple{Vararg{Int}}}: The indices of the subsystems to trace over (can be a single integer or a tuple of integers).
  • dims:::Tuple{Vararg{Int}}: The dimensions of the subsystems.

Returns

The partial trace of the input density matrix. ```

source
QuantumUtilities.partial_traceMethod
partial_trace(v::AbstractVector, trace_indices::Union{Int,Tuple{Vararg{Int}}}, dims::Tuple{Vararg{Int}})

Computes the partial trace of a pure state v.

Arguments

  • v::AbstractVector: The input pure state.
  • trace_indices::Union{Int,Tuple{Vararg{Int}}}: The indices of the subsystems to trace over (can be a single integer or a tuple of integers).
  • dims:Tuple{Vararg{Int}}: The dimensions of the subsystems.

Returns

The partial trace of the input state. ```

source
QuantumUtilities.position_operatorFunction
position_operator(ncutoff::Int, noffset::Int=0)

Create the position operator matrix in the Fock basis for a quantum harmonic oscillator.

Arguments

  • ncutoff::Int: The cutoff value for the number of excitations to be included.
  • noffset::Int=0: The offset for the number of excitations, defaults to 0.

Returns

  • The matrix representation of the position operator.
source
QuantumUtilities.realifcloseMethod
realifclose(x::Complex{T}; tol=eps(T)) where T <: AbstractFloat

Returns the real part of the complex number x if the imaginary part is close to zero within a specified tolerance, and returns x otherwise.

Arguments

  • x::Complex{T}: The input complex number.
  • tol::AbstractFloat: The tolerance for considering the imaginary part close to zero. Defaults to eps(T).

Returns

The real part of x if the imaginary part is close to zero within the specified tolerance, otherwise x itself.

Examples

julia> realifclose(2 + 0im)
2

julia> realifclose(1e-21 + 1e-21im)
1.0e-21

julia> realifclose(1e-21 + 1e-21im; tol=1e-22)
1.0e-21 + 1.0e-21im
source
QuantumUtilities.right_superopMethod
right_superop(A)
right_superop(A, d::Int)

Computes the Liouville space right superoperator representation of a matrix A.

Arguments

  • A: The input matrix.
  • d: (Optional) The size of the identity matrix. If not provided, it is set to the number of rows in A.

Returns

The right superoperator matrix obtained by performing a Kronecker product between the transpose of A and the identity matrix.

Examples

julia> A = [1 2; 3 4]
2×2 Matrix{Int64}:
 1  2
 3  4

julia> right_superop(A)
4×4 Matrix{Int64}:
 1  0  3  0
 0  1  0  3
 2  0  4  0
 0  2  0  4
source
QuantumUtilities.rotation_operatorMethod
rotation_operator(S0::SpinLength, θ, φ, α)

Create the rotation operator matrix for a given spin length and rotation parameters.

Arguments

  • S0::SpinLength: The spin length.
  • θ: The polar angle characterising the direction of the rotation.
  • φ: The azimuthal angle characterising the direction of the rotation.
  • α: The rotation angle.

Returns

  • The matrix representation of the rotation operator for spin S0.
source
QuantumUtilities.rotation_operatorMethod
rotation_operator(S0::SpinLength, n, α)

Create the rotation operator matrix for a given spin length and rotation parameters.

Arguments

  • S0::SpinLength: The spin length.
  • n: The axis of rotation (must be a three-component vector-like).
  • α: The rotation angle.

Returns

  • The matrix representation of the rotation operator for spin S0.
source
QuantumUtilities.s2_operatorMethod
s2_operator(S0::SpinLength)

Create the matrix representation of the squared spin operator in the standard z basis.

Arguments

  • S0::SpinLength: The spin length.

Returns

  • The matrix representation of the $S^2$ operator.
source
QuantumUtilities.scrapMethod
scrap(x::AbstractFloat; tol=eps(typeof(x)))

Returns a floating-point number with a small value set to zero within a specified tolerance.

Arguments

  • x::AbstractFloat: The input floating-point number.
  • tol: The tolerance for considering the input value close to zero. Defaults to eps(typeof(x)).

Returns

A floating-point number with a small value set to zero within the specified tolerance.

Examples

julia> scrap(1.2e-21)
0.0

julia> scrap(1.2e-21; tol=1e-22)
1.2e-21
source
QuantumUtilities.scrapMethod
scrap(x::Complex{T}; tol=eps(T)) where T <: AbstractFloat

Returns a complex number with small real and imaginary parts set to zero within a specified tolerance.

Arguments

  • x::Complex{T}: The input complex number.
  • tol::AbstractFloat: The tolerance for considering the real and imaginary parts close to zero. Defaults to eps(T).

Returns

A complex number with small real and imaginary parts set to zero within the specified tolerance.

Examples

julia> scrap(1e-20 + 1e-21im)
1.0e-20 + 1.0e-21im

julia> scrap(1e-20 + 1e-25im)
1.0e-20
source
QuantumUtilities.sm_operatorMethod
sm_operator(S0::SpinLength)

Create the matrix representation of the lowering spin ladder operator in the standard z basis.

Arguments

  • S0::SpinLength: The spin length.

Returns

  • The matrix representation of the $S_-$ operator.
source
QuantumUtilities.sp_operatorMethod
sp_operator(S0::SpinLength)

Create the matrix representation of the raising spin ladder operator in the standard z basis.

Arguments

  • S0::SpinLength: The spin length.

Returns

  • The matrix representation of the $S_+$ operator.
source
QuantumUtilities.spin_projectionsMethod
spin_projections(S0::SpinLength; rev=false)

Return an iterator over the projections of a spin of length S0. By default, the iterator ranges from smallest to largest projection. The optional argument rev can be used to reverse the order (i.e. from largest to smallest).

Arguments

  • S0::SpinLength: The spin length.
  • rev=false: Whether to reverse the order of the projections.

Returns

  • The iterator over the spin projections.
source
QuantumUtilities.sx_operatorMethod
sx_operator(S0::SpinLength)

Create the matrix representation of the x spin component operator in the standard z basis.

Arguments

  • S0::SpinLength: The spin length.

Returns

  • The matrix representation of the $S_x$ operator.
source
QuantumUtilities.sy_operatorMethod
sy_operator(S0::SpinLength)

Create the matrix representation of the y spin component operator in the standard z basis.

Arguments

  • S0::SpinLength: The spin length.

Returns

  • The matrix representation of the $S_y$ operator.
source
QuantumUtilities.sz_operatorMethod
sz_operator(S0::SpinLength)

Create the matrix representation of the z spin component operator in the standard z basis.

Arguments

  • S0::SpinLength: The spin length.

Returns

  • The matrix representation of the $S_z$ operator.
source
QuantumUtilities.tensorFunction
tensor(A, B...)

Computes the tensor product of matrices or array-like objects A and B.

Arguments

  • A: The first matrix or array-like object.
  • B...: Additional matrices or array-like objects to compute their tensor product with A.

Returns

The tensor product of the input matrices or array-like objects. ```

source
QuantumUtilities.time_evolution_superopMethod
time_evolution_superop(H, t)

Computes the superoperator corresponding to the Hamiltonian evolution of a system governed by the Hamiltonian matrix H over a time t.

Arguments

  • H: The Hamiltonian matrix representing the dynamics of the system.
  • t: The time step for the evolution.

Returns

The superoperator matrix representing the Hamiltonian evolution over the specified time step. ```

source
QuantumUtilities.usincMethod
usinc(x)

Computes the unnormalized sinc function value for the input x, defined as sin(x)/x.

Arguments

  • x: The input value.

Returns

The unnormalized sinc function value of x.

Examples

julia> usinc(0.5)
0.958851077208406

julia> usinc(0.0)
1.0
source
QuantumUtilities.vector_to_operatorMethod
vector_to_operator(v)
vector_to_operator(v, d::Int)

Reshapes a one-dimensional vector v into a matrix of size d by d.

Arguments

  • v: The input one-dimensional vector.
  • d: (Optional) The desired size of the resulting matrix. If not provided, it is calculated as the rounded integer square root of the number of elements in v.

Returns

A matrix of size d by d representing the elements of v reshaped accordingly.

Examples

julia> v = [1, 3, 2, 4]
4-element Vector{Int64}:
 1
 3
 2
 4

julia> vector_to_operator(v)
2×2 Matrix{Int64}:
 1  2
 3  4
source