API
Contents
Quantum object (Qobj) and type
QuantumToolbox.BraQuantumObject
— TypeBraQuantumObject <: QuantumObjectType
Constructor representing a bra state $\langle\psi|$.
QuantumToolbox.Bra
— Constantconst Bra = BraQuantumObject()
A constant representing the type of BraQuantumObject
: a bra state $\langle\psi|$
QuantumToolbox.KetQuantumObject
— TypeKetQuantumObject <: QuantumObjectType
Constructor representing a ket state $|\psi\rangle$.
QuantumToolbox.Ket
— Constantconst Ket = KetQuantumObject()
A constant representing the type of KetQuantumObject
: a ket state $|\psi\rangle$
QuantumToolbox.OperatorQuantumObject
— TypeOperatorQuantumObject <: QuantumObjectType
Constructor representing an operator $\hat{O}$.
QuantumToolbox.Operator
— Constantconst Operator = OperatorQuantumObject()
A constant representing the type of OperatorQuantumObject
: an operator $\hat{O}$
QuantumToolbox.OperatorBraQuantumObject
— TypeOperatorBraQuantumObject <: QuantumObjectType
Constructor representing a bra state in the SuperOperator
formalism $\langle\langle\rho|$.
QuantumToolbox.OperatorBra
— Constantconst OperatorBra = OperatorBraQuantumObject()
A constant representing the type of OperatorBraQuantumObject
: a bra state in the SuperOperator
formalism $\langle\langle\rho|$.
QuantumToolbox.OperatorKetQuantumObject
— TypeOperatorKetQuantumObject <: QuantumObjectType
Constructor representing a ket state in the SuperOperator
formalism $|\rho\rangle\rangle$.
QuantumToolbox.OperatorKet
— Constantconst OperatorKet = OperatorKetQuantumObject()
A constant representing the type of OperatorKetQuantumObject
: a ket state in the SuperOperator
formalism $|\rho\rangle\rangle$
QuantumToolbox.SuperOperatorQuantumObject
— TypeSuperOperatorQuantumObject <: QuantumObjectType
Constructor representing a super-operator $\hat{\mathcal{O}}$ acting on vectorized density operator matrices.
QuantumToolbox.SuperOperator
— Constantconst SuperOperator = SuperOperatorQuantumObject()
A constant representing the type of SuperOperatorQuantumObject
: a super-operator $\hat{\mathcal{O}}$ acting on vectorized density operator matrices
QuantumToolbox.QuantumObject
— Typestruct QuantumObject{MT<:AbstractArray,ObjType<:QuantumObjectType,N}
+API · QuantumToolbox.jl API
Contents
Quantum object (Qobj) and type
QuantumToolbox.BraQuantumObject
— TypeBraQuantumObject <: QuantumObjectType
Constructor representing a bra state $\langle\psi|$.
sourceQuantumToolbox.Bra
— Constantconst Bra = BraQuantumObject()
A constant representing the type of BraQuantumObject
: a bra state $\langle\psi|$
sourceQuantumToolbox.KetQuantumObject
— TypeKetQuantumObject <: QuantumObjectType
Constructor representing a ket state $|\psi\rangle$.
sourceQuantumToolbox.Ket
— Constantconst Ket = KetQuantumObject()
A constant representing the type of KetQuantumObject
: a ket state $|\psi\rangle$
sourceQuantumToolbox.OperatorQuantumObject
— TypeOperatorQuantumObject <: QuantumObjectType
Constructor representing an operator $\hat{O}$.
sourceQuantumToolbox.Operator
— Constantconst Operator = OperatorQuantumObject()
A constant representing the type of OperatorQuantumObject
: an operator $\hat{O}$
sourceQuantumToolbox.OperatorBraQuantumObject
— TypeOperatorBraQuantumObject <: QuantumObjectType
Constructor representing a bra state in the SuperOperator
formalism $\langle\langle\rho|$.
sourceQuantumToolbox.OperatorBra
— Constantconst OperatorBra = OperatorBraQuantumObject()
A constant representing the type of OperatorBraQuantumObject
: a bra state in the SuperOperator
formalism $\langle\langle\rho|$.
sourceQuantumToolbox.OperatorKetQuantumObject
— TypeOperatorKetQuantumObject <: QuantumObjectType
Constructor representing a ket state in the SuperOperator
formalism $|\rho\rangle\rangle$.
sourceQuantumToolbox.OperatorKet
— Constantconst OperatorKet = OperatorKetQuantumObject()
A constant representing the type of OperatorKetQuantumObject
: a ket state in the SuperOperator
formalism $|\rho\rangle\rangle$
sourceQuantumToolbox.SuperOperatorQuantumObject
— TypeSuperOperatorQuantumObject <: QuantumObjectType
Constructor representing a super-operator $\hat{\mathcal{O}}$ acting on vectorized density operator matrices.
sourceQuantumToolbox.SuperOperator
— Constantconst SuperOperator = SuperOperatorQuantumObject()
A constant representing the type of SuperOperatorQuantumObject
: a super-operator $\hat{\mathcal{O}}$ acting on vectorized density operator matrices
sourceQuantumToolbox.QuantumObject
— Typestruct QuantumObject{MT<:AbstractArray,ObjType<:QuantumObjectType,N}
data::MT
type::ObjType
dims::SVector{N, Int}
@@ -13,10 +13,10 @@
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢⎦
julia> a isa QuantumObject
-true
sourceQuantumToolbox.OperatorSum
— Typestruct OperatorSum
A constructor to represent a sum of operators $\sum_i c_i \hat{O}_i$ with a list of coefficients $c_i$ and a list of operators $\hat{O}_i$.
This is very useful when we have to update only the coefficients, without allocating memory by performing the sum of the operators.
sourceBase.size
— Functionsize(A::QuantumObject)
-size(A::QuantumObject, idx::Int)
Returns a tuple containing each dimensions of the array in the QuantumObject
.
Optionally, you can specify an index (idx
) to just get the corresponding dimension of the array.
sourceBase.eltype
— Functioneltype(A::QuantumObject)
Returns the elements type of the matrix or vector corresponding to the QuantumObject
A
.
sourceBase.length
— Functionlength(A::QuantumObject)
Returns the length of the matrix or vector corresponding to the QuantumObject
A
.
sourceQuantumToolbox.isbra
— Functionisbra(A::QuantumObject)
Checks if the QuantumObject
A
is a BraQuantumObject
.
sourceQuantumToolbox.isket
— Functionisket(A::QuantumObject)
Checks if the QuantumObject
A
is a KetQuantumObject
.
sourceQuantumToolbox.isoper
— Functionisoper(A::QuantumObject)
Checks if the QuantumObject
A
is a OperatorQuantumObject
.
sourceQuantumToolbox.isoperbra
— Functionisoperbra(A::QuantumObject)
Checks if the QuantumObject
A
is a OperatorBraQuantumObject
.
sourceQuantumToolbox.isoperket
— Functionisoperket(A::QuantumObject)
Checks if the QuantumObject
A
is a OperatorKetQuantumObject
.
sourceQuantumToolbox.issuper
— Functionissuper(A::QuantumObject)
Checks if the QuantumObject
A
is a SuperOperatorQuantumObject
.
sourceLinearAlgebra.ishermitian
— Functionishermitian(A::QuantumObject)
Test whether the QuantumObject
is Hermitian.
sourceLinearAlgebra.issymmetric
— Functionissymmetric(A::QuantumObject)
Test whether the QuantumObject
is symmetric.
sourceLinearAlgebra.isposdef
— Functionisposdef(A::QuantumObject)
Test whether the QuantumObject
is positive definite (and Hermitian) by trying to perform a Cholesky factorization of A
.
sourceQuantumToolbox.isunitary
— Functionisunitary(U::QuantumObject; kwargs...)
Test whether the QuantumObject
$U$ is unitary operator. This function calls Base.isapprox
to test whether $U U^\dagger$ is approximately equal to identity operator.
Note that all the keyword arguments will be passed to Base.isapprox
.
sourceQobj arithmetic and attributes
Base.conj
— Functionconj(A::QuantumObject)
Return the element-wise complex conjugation of the QuantumObject
.
sourceBase.transpose
— Functiontranspose(A::QuantumObject)
Lazy matrix transpose of the QuantumObject
.
sourceBase.adjoint
— FunctionA'
-adjoint(A::QuantumObject)
Lazy adjoint (conjugate transposition) of the QuantumObject
Note that A'
is a synonym for adjoint(A)
sourceLinearAlgebra.dot
— Functiondot(A::QuantumObject, B::QuantumObject)
Compute the dot product between two QuantumObject
: $\langle A | B \rangle$
Note that A
and B
should be Ket
or OperatorKet
A ⋅ B
(where ⋅
can be typed by tab-completing \cdot
in the REPL) is a synonym for dot(A, B)
sourcedot(i::QuantumObject, A::QuantumObject j::QuantumObject)
Compute the generalized dot product dot(i, A*j)
between three QuantumObject
: $\langle i | \hat{A} | j \rangle$
Supports the following inputs:
A
is in the type of Operator
, with i
and j
are both Ket
.A
is in the type of SuperOperator
, with i
and j
are both OperatorKet
sourceBase.sqrt
— Function√(A)
-sqrt(A::QuantumObject)
Matrix square root of QuantumObject
Note that √(A)
is a synonym for sqrt(A)
sourceBase.log
— Functionlog(A::QuantumObject)
Matrix logarithm of QuantumObject
Note that this function only supports for Operator
and SuperOperator
sourceBase.exp
— Functionexp(A::QuantumObject)
Matrix exponential of QuantumObject
Note that this function only supports for Operator
and SuperOperator
sourceBase.sin
— Functionsin(A::QuantumObject)
Matrix sine of QuantumObject
, defined as
$\sin \left( \hat{A} \right) = \frac{e^{i \hat{A}} - e^{-i \hat{A}}}{2 i}$
Note that this function only supports for Operator
and SuperOperator
sourceBase.cos
— Functioncos(A::QuantumObject)
Matrix cosine of QuantumObject
, defined as
$\cos \left( \hat{A} \right) = \frac{e^{i \hat{A}} + e^{-i \hat{A}}}{2}$
Note that this function only supports for Operator
and SuperOperator
sourceLinearAlgebra.tr
— Functiontr(A::QuantumObject)
Returns the trace of QuantumObject
.
Note that this function only supports for Operator
and SuperOperator
Examples
julia> a = destroy(20)
+true
sourceQuantumToolbox.OperatorSum
— Typestruct OperatorSum
A constructor to represent a sum of operators $\sum_i c_i \hat{O}_i$ with a list of coefficients $c_i$ and a list of operators $\hat{O}_i$.
This is very useful when we have to update only the coefficients, without allocating memory by performing the sum of the operators.
sourceBase.size
— Functionsize(A::QuantumObject)
+size(A::QuantumObject, idx::Int)
Returns a tuple containing each dimensions of the array in the QuantumObject
.
Optionally, you can specify an index (idx
) to just get the corresponding dimension of the array.
sourceBase.eltype
— Functioneltype(A::QuantumObject)
Returns the elements type of the matrix or vector corresponding to the QuantumObject
A
.
sourceBase.length
— Functionlength(A::QuantumObject)
Returns the length of the matrix or vector corresponding to the QuantumObject
A
.
sourceQuantumToolbox.isbra
— Functionisbra(A::QuantumObject)
Checks if the QuantumObject
A
is a BraQuantumObject
.
sourceQuantumToolbox.isket
— Functionisket(A::QuantumObject)
Checks if the QuantumObject
A
is a KetQuantumObject
.
sourceQuantumToolbox.isoper
— Functionisoper(A::QuantumObject)
Checks if the QuantumObject
A
is a OperatorQuantumObject
.
sourceQuantumToolbox.isoperbra
— Functionisoperbra(A::QuantumObject)
Checks if the QuantumObject
A
is a OperatorBraQuantumObject
.
sourceQuantumToolbox.isoperket
— Functionisoperket(A::QuantumObject)
Checks if the QuantumObject
A
is a OperatorKetQuantumObject
.
sourceQuantumToolbox.issuper
— Functionissuper(A::QuantumObject)
Checks if the QuantumObject
A
is a SuperOperatorQuantumObject
.
sourceLinearAlgebra.ishermitian
— Functionishermitian(A::QuantumObject)
Test whether the QuantumObject
is Hermitian.
sourceLinearAlgebra.issymmetric
— Functionissymmetric(A::QuantumObject)
Test whether the QuantumObject
is symmetric.
sourceLinearAlgebra.isposdef
— Functionisposdef(A::QuantumObject)
Test whether the QuantumObject
is positive definite (and Hermitian) by trying to perform a Cholesky factorization of A
.
sourceQuantumToolbox.isunitary
— Functionisunitary(U::QuantumObject; kwargs...)
Test whether the QuantumObject
$U$ is unitary operator. This function calls Base.isapprox
to test whether $U U^\dagger$ is approximately equal to identity operator.
Note that all the keyword arguments will be passed to Base.isapprox
.
sourceQobj arithmetic and attributes
Base.conj
— Functionconj(A::QuantumObject)
Return the element-wise complex conjugation of the QuantumObject
.
sourceBase.transpose
— Functiontranspose(A::QuantumObject)
Lazy matrix transpose of the QuantumObject
.
sourceBase.adjoint
— FunctionA'
+adjoint(A::QuantumObject)
Lazy adjoint (conjugate transposition) of the QuantumObject
Note that A'
is a synonym for adjoint(A)
sourceLinearAlgebra.dot
— Functiondot(A::QuantumObject, B::QuantumObject)
Compute the dot product between two QuantumObject
: $\langle A | B \rangle$
Note that A
and B
should be Ket
or OperatorKet
A ⋅ B
(where ⋅
can be typed by tab-completing \cdot
in the REPL) is a synonym for dot(A, B)
sourcedot(i::QuantumObject, A::QuantumObject j::QuantumObject)
Compute the generalized dot product dot(i, A*j)
between three QuantumObject
: $\langle i | \hat{A} | j \rangle$
Supports the following inputs:
A
is in the type of Operator
, with i
and j
are both Ket
.A
is in the type of SuperOperator
, with i
and j
are both OperatorKet
sourceBase.sqrt
— Function√(A)
+sqrt(A::QuantumObject)
Matrix square root of QuantumObject
Note that √(A)
is a synonym for sqrt(A)
sourceBase.log
— Functionlog(A::QuantumObject)
Matrix logarithm of QuantumObject
Note that this function only supports for Operator
and SuperOperator
sourceBase.exp
— Functionexp(A::QuantumObject)
Matrix exponential of QuantumObject
Note that this function only supports for Operator
and SuperOperator
sourceBase.sin
— Functionsin(A::QuantumObject)
Matrix sine of QuantumObject
, defined as
$\sin \left( \hat{A} \right) = \frac{e^{i \hat{A}} - e^{-i \hat{A}}}{2 i}$
Note that this function only supports for Operator
and SuperOperator
sourceBase.cos
— Functioncos(A::QuantumObject)
Matrix cosine of QuantumObject
, defined as
$\cos \left( \hat{A} \right) = \frac{e^{i \hat{A}} + e^{-i \hat{A}}}{2}$
Note that this function only supports for Operator
and SuperOperator
sourceLinearAlgebra.tr
— Functiontr(A::QuantumObject)
Returns the trace of QuantumObject
.
Note that this function only supports for Operator
and SuperOperator
Examples
julia> a = destroy(20)
Quantum Object: type=Operator dims=[20] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 19 stored entries:
⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀
@@ -26,7 +26,7 @@
⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢
julia> tr(a' * a)
-190.0 + 0.0im
sourceLinearAlgebra.svdvals
— Functionsvdvals(A::QuantumObject)
Return the singular values of a QuantumObject
in descending order
sourceLinearAlgebra.norm
— Functionnorm(A::QuantumObject, p::Real)
Return the standard vector p
-norm or Schatten p
-norm of a QuantumObject
depending on the type of A
:
- If
A
is either Ket
, Bra
, OperatorKet
, or OperatorBra
, returns the standard vector p
-norm (default p=2
) of A
. - If
A
is either Operator
or SuperOperator
, returns Schatten p
-norm (default p=1
) of A
.
Examples
julia> ψ = fock(10, 2)
+190.0 + 0.0im
sourceLinearAlgebra.svdvals
— Functionsvdvals(A::QuantumObject)
Return the singular values of a QuantumObject
in descending order
sourceLinearAlgebra.norm
— Functionnorm(A::QuantumObject, p::Real)
Return the standard vector p
-norm or Schatten p
-norm of a QuantumObject
depending on the type of A
:
- If
A
is either Ket
, Bra
, OperatorKet
, or OperatorBra
, returns the standard vector p
-norm (default p=2
) of A
. - If
A
is either Operator
or SuperOperator
, returns Schatten p
-norm (default p=1
) of A
.
Examples
julia> ψ = fock(10, 2)
Quantum Object: type=Ket dims=[10] size=(10,)
10-element Vector{ComplexF64}:
0.0 + 0.0im
@@ -41,7 +41,7 @@
0.0 + 0.0im
julia> norm(ψ)
-1.0
sourceLinearAlgebra.normalize
— Functionnormalize(A::QuantumObject, p::Real)
Return normalized QuantumObject
so that its p
-norm equals to unity, i.e. norm(A, p) == 1
.
Support for the following types of QuantumObject
:
Also, see norm
about its definition for different types of QuantumObject
.
sourceLinearAlgebra.normalize!
— Functionnormalize!(A::QuantumObject, p::Real)
Normalize QuantumObject
in-place so that its p
-norm equals to unity, i.e. norm(A, p) == 1
.
Support for the following types of QuantumObject
:
Also, see norm
about its definition for different types of QuantumObject
.
sourceBase.inv
— Functioninv(A::QuantumObject)
Matrix inverse of the QuantumObject
sourceLinearAlgebra.diag
— Functiondiag(A::QuantumObject, k::Int=0)
Return the k
-th diagonal elements of a matrix-type QuantumObject
Note that this function only supports for Operator
and SuperOperator
sourceQuantumToolbox.proj
— Functionproj(ψ::QuantumObject)
Return the projector for a Ket
or Bra
type of QuantumObject
sourceQuantumToolbox.ptrace
— Functionptrace(QO::QuantumObject, sel)
Partial trace of a quantum state QO
leaving only the dimensions with the indices present in the sel
vector.
Note that this function will always return Operator
. No matter the input QuantumObject
is a Ket
, Bra
, or Operator
.
Examples
Two qubits in the state $\ket{\psi} = \ket{e,g}$:
julia> ψ = kron(fock(2,0), fock(2,1))
+1.0
sourceLinearAlgebra.normalize
— Functionnormalize(A::QuantumObject, p::Real)
Return normalized QuantumObject
so that its p
-norm equals to unity, i.e. norm(A, p) == 1
.
Support for the following types of QuantumObject
:
Also, see norm
about its definition for different types of QuantumObject
.
sourceLinearAlgebra.normalize!
— Functionnormalize!(A::QuantumObject, p::Real)
Normalize QuantumObject
in-place so that its p
-norm equals to unity, i.e. norm(A, p) == 1
.
Support for the following types of QuantumObject
:
Also, see norm
about its definition for different types of QuantumObject
.
sourceBase.inv
— Functioninv(A::QuantumObject)
Matrix inverse of the QuantumObject
sourceLinearAlgebra.diag
— Functiondiag(A::QuantumObject, k::Int=0)
Return the k
-th diagonal elements of a matrix-type QuantumObject
Note that this function only supports for Operator
and SuperOperator
sourceQuantumToolbox.proj
— Functionproj(ψ::QuantumObject)
Return the projector for a Ket
or Bra
type of QuantumObject
sourceQuantumToolbox.ptrace
— Functionptrace(QO::QuantumObject, sel)
Partial trace of a quantum state QO
leaving only the dimensions with the indices present in the sel
vector.
Note that this function will always return Operator
. No matter the input QuantumObject
is a Ket
, Bra
, or Operator
.
Examples
Two qubits in the state $\ket{\psi} = \ket{e,g}$:
julia> ψ = kron(fock(2,0), fock(2,1))
Quantum Object: type=Ket dims=[2, 2] size=(4,)
4-element Vector{ComplexF64}:
0.0 + 0.0im
@@ -65,12 +65,12 @@
Quantum Object: type=Operator dims=[2] size=(2, 2) ishermitian=true
2×2 Matrix{ComplexF64}:
0.5+0.0im 0.0+0.0im
- 0.0+0.0im 0.5+0.0im
sourceQuantumToolbox.purity
— Functionpurity(ρ::QuantumObject)
Calculate the purity of a QuantumObject
: $\textrm{Tr}(\rho^2)$
Note that this function only supports for Ket
, Bra
, and Operator
sourceQuantumToolbox.permute
— Functionpermute(A::QuantumObject, order::Union{AbstractVector{Int},Tuple})
Permute the tensor structure of a QuantumObject
A
according to the specified order
list
Note that this method currently works for Ket
, Bra
, and Operator
types of QuantumObject
.
Examples
If order = [2, 1, 3]
, the Hilbert space structure will be re-arranged: $\mathcal{H}_1 \otimes \mathcal{H}_2 \otimes \mathcal{H}_3 \rightarrow \mathcal{H}_2 \otimes \mathcal{H}_1 \otimes \mathcal{H}_3$.
julia> ψ1 = fock(2, 0)
+ 0.0+0.0im 0.5+0.0im
sourceQuantumToolbox.purity
— Functionpurity(ρ::QuantumObject)
Calculate the purity of a QuantumObject
: $\textrm{Tr}(\rho^2)$
Note that this function only supports for Ket
, Bra
, and Operator
sourceQuantumToolbox.permute
— Functionpermute(A::QuantumObject, order::Union{AbstractVector{Int},Tuple})
Permute the tensor structure of a QuantumObject
A
according to the specified order
list
Note that this method currently works for Ket
, Bra
, and Operator
types of QuantumObject
.
Examples
If order = [2, 1, 3]
, the Hilbert space structure will be re-arranged: $\mathcal{H}_1 \otimes \mathcal{H}_2 \otimes \mathcal{H}_3 \rightarrow \mathcal{H}_2 \otimes \mathcal{H}_1 \otimes \mathcal{H}_3$.
julia> ψ1 = fock(2, 0)
julia> ψ2 = fock(3, 1)
julia> ψ3 = fock(4, 2)
julia> ψ_123 = tensor(ψ1, ψ2, ψ3)
julia> permute(ψ_123, [2, 1, 3]) ≈ tensor(ψ2, ψ1, ψ3)
-true
Beware of type-stability! It is highly recommended to use permute(A, order)
with order
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.tidyup
— Functiontidyup(A::QuantumObject, tol::Real=1e-14)
Given a QuantumObject
A
, check the real and imaginary parts of each element separately. Remove the real or imaginary value if its absolute value is less than tol
.
sourceQuantumToolbox.tidyup!
— Functiontidyup!(A::QuantumObject, tol::Real=1e-14)
Given a QuantumObject
A
, check the real and imaginary parts of each element separately. Remove the real or imaginary value if its absolute value is less than tol
.
Note that this function is an in-place version of tidyup
.
sourceQuantumToolbox.get_data
— Functionget_data(A::QuantumObject)
Returns the data of a QuantumObject.
sourceQuantumToolbox.get_coherence
— Functionget_coherence(ψ::QuantumObject)
Get the coherence value $\alpha$ by measuring the expectation value of the destruction operator $\hat{a}$ on a state $\ket{\psi}$ or a density matrix $\hat{\rho}$.
It returns both $\alpha$ and the corresponding state with the coherence removed: $\ket{\delta_\alpha} = \exp ( \alpha^* \hat{a} - \alpha \hat{a}^\dagger ) \ket{\psi}$ for a pure state, and $\hat{\rho_\alpha} = \exp ( \alpha^* \hat{a} - \alpha \hat{a}^\dagger ) \hat{\rho} \exp ( -\bar{\alpha} \hat{a} + \alpha \hat{a}^\dagger )$ for a density matrix. These states correspond to the quantum fluctuations around the coherent state $\ket{\alpha}$ or $|\alpha\rangle\langle\alpha|$.
sourceQuantumToolbox.partial_transpose
— Functionpartial_transpose(ρ::QuantumObject, mask::Vector{Bool})
Return the partial transpose of a density matrix $\rho$, where mask
is an array/vector with length that equals the length of ρ.dims
. The elements in mask
are boolean (true
or false
) which indicates whether or not the corresponding subsystem should be transposed.
Arguments
ρ::QuantumObject
: The density matrix (ρ.type
must be OperatorQuantumObject
).mask::Vector{Bool}
: A boolean vector selects which subsystems should be transposed.
Returns
ρ_pt::QuantumObject
: The density matrix with the selected subsystems transposed.
sourceQobj eigenvalues and eigenvectors
QuantumToolbox.EigsolveResult
— Typestruct EigsolveResult{T1<:Vector{<:Number}, T2<:AbstractMatrix{<:Number}, ObjType<:Union{Nothing,OperatorQuantumObject,SuperOperatorQuantumObject},N}
+true
Beware of type-stability! It is highly recommended to use permute(A, order)
with order
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.tidyup
— Functiontidyup(A::QuantumObject, tol::Real=1e-14)
Given a QuantumObject
A
, check the real and imaginary parts of each element separately. Remove the real or imaginary value if its absolute value is less than tol
.
sourceQuantumToolbox.tidyup!
— Functiontidyup!(A::QuantumObject, tol::Real=1e-14)
Given a QuantumObject
A
, check the real and imaginary parts of each element separately. Remove the real or imaginary value if its absolute value is less than tol
.
Note that this function is an in-place version of tidyup
.
sourceQuantumToolbox.get_data
— Functionget_data(A::QuantumObject)
Returns the data of a QuantumObject.
sourceQuantumToolbox.get_coherence
— Functionget_coherence(ψ::QuantumObject)
Get the coherence value $\alpha$ by measuring the expectation value of the destruction operator $\hat{a}$ on a state $\ket{\psi}$ or a density matrix $\hat{\rho}$.
It returns both $\alpha$ and the corresponding state with the coherence removed: $\ket{\delta_\alpha} = \exp ( \alpha^* \hat{a} - \alpha \hat{a}^\dagger ) \ket{\psi}$ for a pure state, and $\hat{\rho_\alpha} = \exp ( \alpha^* \hat{a} - \alpha \hat{a}^\dagger ) \hat{\rho} \exp ( -\bar{\alpha} \hat{a} + \alpha \hat{a}^\dagger )$ for a density matrix. These states correspond to the quantum fluctuations around the coherent state $\ket{\alpha}$ or $|\alpha\rangle\langle\alpha|$.
sourceQuantumToolbox.partial_transpose
— Functionpartial_transpose(ρ::QuantumObject, mask::Vector{Bool})
Return the partial transpose of a density matrix $\rho$, where mask
is an array/vector with length that equals the length of ρ.dims
. The elements in mask
are boolean (true
or false
) which indicates whether or not the corresponding subsystem should be transposed.
Arguments
ρ::QuantumObject
: The density matrix (ρ.type
must be OperatorQuantumObject
).mask::Vector{Bool}
: A boolean vector selects which subsystems should be transposed.
Returns
ρ_pt::QuantumObject
: The density matrix with the selected subsystems transposed.
sourceQobj eigenvalues and eigenvectors
QuantumToolbox.EigsolveResult
— Typestruct EigsolveResult{T1<:Vector{<:Number}, T2<:AbstractMatrix{<:Number}, ObjType<:Union{Nothing,OperatorQuantumObject,SuperOperatorQuantumObject},N}
values::T1
vectors::T2
type::ObjType
@@ -95,7 +95,7 @@
julia> T
2×2 Matrix{ComplexF64}:
-0.707107+0.0im 0.707107+0.0im
- 0.707107+0.0im 0.707107+0.0im
sourceQuantumToolbox.eigenenergies
— Functioneigenenergies(A::QuantumObject; sparse::Bool=false, kwargs...)
Calculate the eigenenergies
Arguments
A::QuantumObject
: the QuantumObject
to solve eigenvaluessparse::Bool
: if false
call eigvals(A::QuantumObject; kwargs...)
, otherwise call eigsolve
. Default to false
.kwargs
: Additional keyword arguments passed to the solver
Returns
::Vector{<:Number}
: a list of eigenvalues
sourceQuantumToolbox.eigenstates
— Functioneigenstates(A::QuantumObject; sparse::Bool=false, kwargs...)
Calculate the eigenvalues and corresponding eigenvectors
Arguments
A::QuantumObject
: the QuantumObject
to solve eigenvalues and eigenvectorssparse::Bool
: if false
call eigen(A::QuantumObject; kwargs...)
, otherwise call eigsolve
. Default to false
.kwargs
: Additional keyword arguments passed to the solver
Returns
::EigsolveResult
: containing the eigenvalues, the eigenvectors, and some information from the solver. see also EigsolveResult
sourceLinearAlgebra.eigen
— FunctionLinearAlgebra.eigen(A::QuantumObject; kwargs...)
Calculates the eigenvalues and eigenvectors of the QuantumObject
A
using the Julia LinearAlgebra package.
julia> a = destroy(5);
+ 0.707107+0.0im 0.707107+0.0im
sourceQuantumToolbox.eigenenergies
— Functioneigenenergies(A::QuantumObject; sparse::Bool=false, kwargs...)
Calculate the eigenenergies
Arguments
A::QuantumObject
: the QuantumObject
to solve eigenvaluessparse::Bool
: if false
call eigvals(A::QuantumObject; kwargs...)
, otherwise call eigsolve
. Default to false
.kwargs
: Additional keyword arguments passed to the solver
Returns
::Vector{<:Number}
: a list of eigenvalues
sourceQuantumToolbox.eigenstates
— Functioneigenstates(A::QuantumObject; sparse::Bool=false, kwargs...)
Calculate the eigenvalues and corresponding eigenvectors
Arguments
A::QuantumObject
: the QuantumObject
to solve eigenvalues and eigenvectorssparse::Bool
: if false
call eigen(A::QuantumObject; kwargs...)
, otherwise call eigsolve
. Default to false
.kwargs
: Additional keyword arguments passed to the solver
Returns
::EigsolveResult
: containing the eigenvalues, the eigenvectors, and some information from the solver. see also EigsolveResult
sourceLinearAlgebra.eigen
— FunctionLinearAlgebra.eigen(A::QuantumObject; kwargs...)
Calculates the eigenvalues and eigenvectors of the QuantumObject
A
using the Julia LinearAlgebra package.
julia> a = destroy(5);
julia> H = a + a'
Quantum Object: type=Operator dims=[5] size=(5, 5) ishermitian=true
@@ -124,7 +124,7 @@
0.447214+0.0im 0.447214+0.0im -0.447214-0.0im 0.447214-0.0im
julia> expect(H, ψ[1]) ≈ E[1]
-true
sourceLinearAlgebra.eigvals
— FunctionLinearAlgebra.eigvals(A::QuantumObject; kwargs...)
Same as eigen(A::QuantumObject; kwargs...)
but for only the eigenvalues.
sourceQuantumToolbox.eigsolve
— Functioneigsolve(A::QuantumObject;
+true
sourceLinearAlgebra.eigvals
— FunctionLinearAlgebra.eigvals(A::QuantumObject; kwargs...)
Same as eigen(A::QuantumObject; kwargs...)
but for only the eigenvalues.
sourceQuantumToolbox.eigsolve
— Functioneigsolve(A::QuantumObject;
v0::Union{Nothing,AbstractVector}=nothing,
sigma::Union{Nothing, Real}=nothing,
k::Int = 1,
@@ -132,7 +132,7 @@
tol::Real = 1e-8,
maxiter::Int = 200,
solver::Union{Nothing, SciMLLinearSolveAlgorithm} = nothing,
- kwargs...)
Solve for the eigenvalues and eigenvectors of a matrix A
using the Arnoldi method.
Notes
- For more details about
solver
and extra kwargs
, please refer to LinearSolve.jl
Returns
EigsolveResult
: A struct containing the eigenvalues, the eigenvectors, and some information about the eigsolver
sourceQuantumToolbox.eigsolve_al
— Functioneigsolve_al(H::QuantumObject,
+ kwargs...)
Solve for the eigenvalues and eigenvectors of a matrix A
using the Arnoldi method.
Notes
- For more details about
solver
and extra kwargs
, please refer to LinearSolve.jl
Returns
EigsolveResult
: A struct containing the eigenvalues, the eigenvectors, and some information about the eigsolver
sourceQuantumToolbox.eigsolve_al
— Functioneigsolve_al(H::QuantumObject,
T::Real, c_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
alg::OrdinaryDiffEqAlgorithm=Tsit5(),
H_t::Union{Nothing,Function}=nothing,
@@ -142,7 +142,7 @@
krylovdim::Int=min(10, size(H, 1)),
maxiter::Int=200,
eigstol::Real=1e-6,
- kwargs...)
Solve the eigenvalue problem for a Liouvillian superoperator L
using the Arnoldi-Lindblad method.
Arguments
H
: The Hamiltonian (or directly the Liouvillian) of the system.T
: The time at which to evaluate the time evolutionc_ops
: A vector of collapse operators. Default is nothing
meaning the system is closed.alg
: The differential equation solver algorithmH_t
: A function H_t(t)
that returns the additional term at time t
params
: A dictionary of additional parametersρ0
: The initial density matrix. If not specified, a random density matrix is usedk
: The number of eigenvalues to computekrylovdim
: The dimension of the Krylov subspacemaxiter
: The maximum number of iterations for the eigsolvereigstol
: The tolerance for the eigsolverkwargs
: Additional keyword arguments passed to the differential equation solver
Notes
- For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
EigsolveResult
: A struct containing the eigenvalues, the eigenvectors, and some information about the eigsolver
References
- [1] Minganti, F., & Huybrechts, D. (2022). Arnoldi-Lindblad time evolution: Faster-than-the-clock algorithm for the spectrum of time-independent and Floquet open quantum systems. Quantum, 6, 649.
sourceQobj manipulation
QuantumToolbox.ket2dm
— Functionket2dm(ψ::QuantumObject)
Transform the ket state $\ket{\psi}$ into a pure density matrix $\hat{\rho} = \dyad{\psi}$.
sourceQuantumToolbox.expect
— Functionexpect(O::QuantumObject, ψ::Union{QuantumObject,Vector{QuantumObject}})
Expectation value of the Operator
O
with the state ψ
. The state can be a Ket
, Bra
or Operator
.
If ψ
is a Ket
or Bra
, the function calculates $\langle\psi|\hat{O}|\psi\rangle$.
If ψ
is a density matrix (Operator
), the function calculates $\textrm{Tr} \left[ \hat{O} \hat{\psi} \right]$
The function returns a real number if O
is of Hermitian
type or Symmetric
type, and returns a complex number otherwise. You can make an operator O
hermitian by using Hermitian(O)
.
Note that ψ
can also be given as a list of QuantumObject
, it returns a list of expectation values.
Examples
julia> ψ = 1 / √2 * (fock(10,2) + fock(10,4));
+ kwargs...)
Solve the eigenvalue problem for a Liouvillian superoperator L
using the Arnoldi-Lindblad method.
Arguments
H
: The Hamiltonian (or directly the Liouvillian) of the system.T
: The time at which to evaluate the time evolutionc_ops
: A vector of collapse operators. Default is nothing
meaning the system is closed.alg
: The differential equation solver algorithmH_t
: A function H_t(t)
that returns the additional term at time t
params
: A dictionary of additional parametersρ0
: The initial density matrix. If not specified, a random density matrix is usedk
: The number of eigenvalues to computekrylovdim
: The dimension of the Krylov subspacemaxiter
: The maximum number of iterations for the eigsolvereigstol
: The tolerance for the eigsolverkwargs
: Additional keyword arguments passed to the differential equation solver
Notes
- For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
EigsolveResult
: A struct containing the eigenvalues, the eigenvectors, and some information about the eigsolver
References
- [1] Minganti, F., & Huybrechts, D. (2022). Arnoldi-Lindblad time evolution: Faster-than-the-clock algorithm for the spectrum of time-independent and Floquet open quantum systems. Quantum, 6, 649.
sourceQobj manipulation
QuantumToolbox.ket2dm
— Functionket2dm(ψ::QuantumObject)
Transform the ket state $\ket{\psi}$ into a pure density matrix $\hat{\rho} = \dyad{\psi}$.
sourceQuantumToolbox.expect
— Functionexpect(O::QuantumObject, ψ::Union{QuantumObject,Vector{QuantumObject}})
Expectation value of the Operator
O
with the state ψ
. The state can be a Ket
, Bra
or Operator
.
If ψ
is a Ket
or Bra
, the function calculates $\langle\psi|\hat{O}|\psi\rangle$.
If ψ
is a density matrix (Operator
), the function calculates $\textrm{Tr} \left[ \hat{O} \hat{\psi} \right]$
The function returns a real number if O
is of Hermitian
type or Symmetric
type, and returns a complex number otherwise. You can make an operator O
hermitian by using Hermitian(O)
.
Note that ψ
can also be given as a list of QuantumObject
, it returns a list of expectation values.
Examples
julia> ψ = 1 / √2 * (fock(10,2) + fock(10,4));
julia> a = destroy(10);
@@ -150,7 +150,7 @@
3.0 + 0.0im
julia> expect(Hermitian(a' * a), ψ) |> round
-3.0
sourceQuantumToolbox.variance
— Functionvariance(O::QuantumObject, ψ::Union{QuantumObject,Vector{QuantumObject}})
Variance of the Operator
O
: $\langle\hat{O}^2\rangle - \langle\hat{O}\rangle^2$,
where $\langle\hat{O}\rangle$ is the expectation value of O
with the state ψ
(see also expect
), and the state ψ
can be a Ket
, Bra
or Operator
.
The function returns a real number if O
is hermitian, and returns a complex number otherwise.
Note that ψ
can also be given as a list of QuantumObject
, it returns a list of expectation values.
sourceBase.kron
— Functionkron(A::QuantumObject, B::QuantumObject, ...)
Returns the Kronecker product $\hat{A} \otimes \hat{B} \otimes \cdots$.
Examples
julia> a = destroy(20)
+3.0
sourceQuantumToolbox.variance
— Functionvariance(O::QuantumObject, ψ::Union{QuantumObject,Vector{QuantumObject}})
Variance of the Operator
O
: $\langle\hat{O}^2\rangle - \langle\hat{O}\rangle^2$,
where $\langle\hat{O}\rangle$ is the expectation value of O
with the state ψ
(see also expect
), and the state ψ
can be a Ket
, Bra
or Operator
.
The function returns a real number if O
is hermitian, and returns a complex number otherwise.
Note that ψ
can also be given as a list of QuantumObject
, it returns a list of expectation values.
sourceBase.kron
— Functionkron(A::QuantumObject, B::QuantumObject, ...)
Returns the Kronecker product $\hat{A} \otimes \hat{B} \otimes \cdots$.
Examples
julia> a = destroy(20)
Quantum Object: type=Operator dims=[20] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 19 stored entries:
⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀
@@ -182,7 +182,7 @@
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡀⠀
-⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⠦
sourceQuantumToolbox.sparse_to_dense
— Functionsparse_to_dense(A::QuantumObject)
Converts a sparse QuantumObject to a dense QuantumObject.
sourceQuantumToolbox.dense_to_sparse
— Functiondense_to_sparse(A::QuantumObject)
Converts a dense QuantumObject to a sparse QuantumObject.
sourceQuantumToolbox.vec2mat
— Functionvec2mat(A::AbstractVector)
Converts a vector to a matrix.
sourcevec2mat(A::QuantumObject)
Convert a quantum object from vector (OperatorKetQuantumObject
-type) to matrix (OperatorQuantumObject
-type)
sourceQuantumToolbox.mat2vec
— Functionmat2vec(A::QuantumObject)
Convert a quantum object from matrix (OperatorQuantumObject
-type) to vector (OperatorKetQuantumObject
-type)
sourcemat2vec(A::AbstractMatrix)
Converts a matrix to a vector.
sourceGenerate states and operators
QuantumToolbox.zero_ket
— Functionzero_ket(dimensions)
Returns a zero Ket
vector with given argument dimensions
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int}, Tuple}
: list of dimensions representing the each number of basis in the subsystems.
Beware of type-stability! It is highly recommended to use zero_ket(dimensions)
with dimensions
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.fock
— Functionfock(N::Int, j::Int=0; dims::Union{Int,AbstractVector{Int},Tuple}=N, sparse::Union{Bool,Val}=Val(false))
Generates a fock state $\ket{\psi}$ of dimension N
.
It is also possible to specify the list of dimensions dims
if different subsystems are present.
Beware of type-stability! If you want to keep type stability, it is recommended to use fock(N, j, dims=dims, sparse=Val(sparse))
instead of fock(N, j, dims=dims, sparse=sparse)
. Consider also to use dims
as a Tuple
or SVector
instead of Vector
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.basis
— Functionbasis(N::Int, j::Int = 0; dims::Union{Int,AbstractVector{Int},Tuple}=N)
Generates a fock state like fock
.
It is also possible to specify the list of dimensions dims
if different subsystems are present.
Beware of type-stability! If you want to keep type stability, it is recommended to use basis(N, j, dims=dims)
with dims
as a Tuple
or SVector
instead of Vector
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.coherent
— Functioncoherent(N::Int, α::Number)
Generates a coherent state $|\alpha\rangle$, which is defined as an eigenvector of the bosonic annihilation operator $\hat{a} |\alpha\rangle = \alpha |\alpha\rangle$.
This state is constructed via the displacement operator displace
and zero-fock state fock
: $|\alpha\rangle = \hat{D}(\alpha) |0\rangle$
sourceQuantumToolbox.rand_ket
— Functionrand_ket(dimensions)
Generate a random normalized Ket
vector with given argument dimensions
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
Beware of type-stability! If you want to keep type stability, it is recommended to use rand_ket(dimensions)
with dimensions
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.fock_dm
— Functionfock_dm(N::Int, j::Int=0; dims::Union{Int,AbstractVector{Int},Tuple}=N, sparse::Union{Bool,Val}=Val(false))
Density matrix representation of a Fock state.
Constructed via outer product of fock
.
Beware of type-stability! If you want to keep type stability, it is recommended to use fock_dm(N, j, dims=dims, sparse=Val(sparse))
instead of fock_dm(N, j, dims=dims, sparse=sparse)
. Consider also to use dims
as a Tuple
or SVector
instead of Vector
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.coherent_dm
— Functioncoherent_dm(N::Int, α::Number)
Density matrix representation of a coherent state.
Constructed via outer product of coherent
.
sourceQuantumToolbox.thermal_dm
— Functionthermal_dm(N::Int, n::Real; sparse::Union{Bool,Val}=Val(false))
Density matrix for a thermal state (generating thermal state probabilities) with the following arguments:
N::Int
: Number of basis states in the Hilbert spacen::Real
: Expectation value for number of particles in the thermal state.sparse::Union{Bool,Val}
: If true
, return a sparse matrix representation.
Beware of type-stability! If you want to keep type stability, it is recommended to use thermal_dm(N, n, sparse=Val(sparse))
instead of thermal_dm(N, n, sparse=sparse)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.maximally_mixed_dm
— Functionmaximally_mixed_dm(dimensions)
Returns the maximally mixed density matrix with given argument dimensions
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
Beware of type-stability! If you want to keep type stability, it is recommended to use maximally_mixed_dm(dimensions)
with dimensions
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.rand_dm
— Functionrand_dm(dimensions; rank::Int=prod(dimensions))
Generate a random density matrix from Ginibre ensemble with given argument dimensions
and rank
, ensuring that it is positive semi-definite and trace equals to 1
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
The default keyword argument rank = prod(dimensions)
(full rank).
Beware of type-stability! If you want to keep type stability, it is recommended to use rand_dm(dimensions; rank=rank)
with dimensions
as Tuple
or SVector
instead of Vector
. See this link and the related Section about type stability for more details.
References
sourceQuantumToolbox.spin_state
— Functionspin_state(j::Real, m::Real)
Generate the spin state: $|j, m\rangle$
The eigenstate of the Spin-j
$\hat{S}_z$ operator with eigenvalue m
, where where j
is the spin quantum number and can be a non-negative integer or half-integer
See also jmat
.
sourceQuantumToolbox.spin_coherent
— Functionspin_coherent(j::Real, θ::Real, ϕ::Real)
Generate the coherent spin state (rotation of the $|j, j\rangle$ state), namely
\[|\theta, \phi \rangle = \hat{R}(\theta, \phi) |j, j\rangle\]
where the rotation operator is defined as
\[\hat{R}(\theta, \phi) = \exp \left( \frac{\theta}{2} (\hat{S}_- e^{i\phi} - \hat{S}_+ e^{-i\phi}) \right)\]
and $\hat{S}_\pm$ are plus and minus Spin-j
operators, respectively.
Arguments
j::Real
: The spin quantum number and can be a non-negative integer or half-integerθ::Real
: rotation angle from z-axisϕ::Real
: rotation angle from x-axis
See also jmat
and spin_state
.
Reference
sourceQuantumToolbox.bell_state
— Functionbell_state(x::Union{Int}, z::Union{Int})
Return the Bell state depending on the arguments (x, z)
:
(0, 0)
: $| \Phi^+ \rangle = ( |00\rangle + |11\rangle ) / \sqrt{2}$(0, 1)
: $| \Phi^- \rangle = ( |00\rangle - |11\rangle ) / \sqrt{2}$(1, 0)
: $| \Psi^+ \rangle = ( |01\rangle + |10\rangle ) / \sqrt{2}$(1, 1)
: $| \Psi^- \rangle = ( |01\rangle - |10\rangle ) / \sqrt{2}$
Here, x = 1
(z = 1
) means applying Pauli-$X$ ( Pauli-$Z$) unitary transformation on $| \Phi^+ \rangle$.
Example
julia> bell_state(0, 0)
+⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⠦
sourceQuantumToolbox.sparse_to_dense
— Functionsparse_to_dense(A::QuantumObject)
Converts a sparse QuantumObject to a dense QuantumObject.
sourceQuantumToolbox.dense_to_sparse
— Functiondense_to_sparse(A::QuantumObject)
Converts a dense QuantumObject to a sparse QuantumObject.
sourceQuantumToolbox.vec2mat
— Functionvec2mat(A::AbstractVector)
Converts a vector to a matrix.
sourcevec2mat(A::QuantumObject)
Convert a quantum object from vector (OperatorKetQuantumObject
-type) to matrix (OperatorQuantumObject
-type)
sourceQuantumToolbox.mat2vec
— Functionmat2vec(A::QuantumObject)
Convert a quantum object from matrix (OperatorQuantumObject
-type) to vector (OperatorKetQuantumObject
-type)
sourcemat2vec(A::AbstractMatrix)
Converts a matrix to a vector.
sourceGenerate states and operators
QuantumToolbox.zero_ket
— Functionzero_ket(dimensions)
Returns a zero Ket
vector with given argument dimensions
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int}, Tuple}
: list of dimensions representing the each number of basis in the subsystems.
Beware of type-stability! It is highly recommended to use zero_ket(dimensions)
with dimensions
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.fock
— Functionfock(N::Int, j::Int=0; dims::Union{Int,AbstractVector{Int},Tuple}=N, sparse::Union{Bool,Val}=Val(false))
Generates a fock state $\ket{\psi}$ of dimension N
.
It is also possible to specify the list of dimensions dims
if different subsystems are present.
Beware of type-stability! If you want to keep type stability, it is recommended to use fock(N, j, dims=dims, sparse=Val(sparse))
instead of fock(N, j, dims=dims, sparse=sparse)
. Consider also to use dims
as a Tuple
or SVector
instead of Vector
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.basis
— Functionbasis(N::Int, j::Int = 0; dims::Union{Int,AbstractVector{Int},Tuple}=N)
Generates a fock state like fock
.
It is also possible to specify the list of dimensions dims
if different subsystems are present.
Beware of type-stability! If you want to keep type stability, it is recommended to use basis(N, j, dims=dims)
with dims
as a Tuple
or SVector
instead of Vector
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.coherent
— Functioncoherent(N::Int, α::Number)
Generates a coherent state $|\alpha\rangle$, which is defined as an eigenvector of the bosonic annihilation operator $\hat{a} |\alpha\rangle = \alpha |\alpha\rangle$.
This state is constructed via the displacement operator displace
and zero-fock state fock
: $|\alpha\rangle = \hat{D}(\alpha) |0\rangle$
sourceQuantumToolbox.rand_ket
— Functionrand_ket(dimensions)
Generate a random normalized Ket
vector with given argument dimensions
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
Beware of type-stability! If you want to keep type stability, it is recommended to use rand_ket(dimensions)
with dimensions
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.fock_dm
— Functionfock_dm(N::Int, j::Int=0; dims::Union{Int,AbstractVector{Int},Tuple}=N, sparse::Union{Bool,Val}=Val(false))
Density matrix representation of a Fock state.
Constructed via outer product of fock
.
Beware of type-stability! If you want to keep type stability, it is recommended to use fock_dm(N, j, dims=dims, sparse=Val(sparse))
instead of fock_dm(N, j, dims=dims, sparse=sparse)
. Consider also to use dims
as a Tuple
or SVector
instead of Vector
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.coherent_dm
— Functioncoherent_dm(N::Int, α::Number)
Density matrix representation of a coherent state.
Constructed via outer product of coherent
.
sourceQuantumToolbox.thermal_dm
— Functionthermal_dm(N::Int, n::Real; sparse::Union{Bool,Val}=Val(false))
Density matrix for a thermal state (generating thermal state probabilities) with the following arguments:
N::Int
: Number of basis states in the Hilbert spacen::Real
: Expectation value for number of particles in the thermal state.sparse::Union{Bool,Val}
: If true
, return a sparse matrix representation.
Beware of type-stability! If you want to keep type stability, it is recommended to use thermal_dm(N, n, sparse=Val(sparse))
instead of thermal_dm(N, n, sparse=sparse)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.maximally_mixed_dm
— Functionmaximally_mixed_dm(dimensions)
Returns the maximally mixed density matrix with given argument dimensions
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
Beware of type-stability! If you want to keep type stability, it is recommended to use maximally_mixed_dm(dimensions)
with dimensions
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.rand_dm
— Functionrand_dm(dimensions; rank::Int=prod(dimensions))
Generate a random density matrix from Ginibre ensemble with given argument dimensions
and rank
, ensuring that it is positive semi-definite and trace equals to 1
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
The default keyword argument rank = prod(dimensions)
(full rank).
Beware of type-stability! If you want to keep type stability, it is recommended to use rand_dm(dimensions; rank=rank)
with dimensions
as Tuple
or SVector
instead of Vector
. See this link and the related Section about type stability for more details.
References
sourceQuantumToolbox.spin_state
— Functionspin_state(j::Real, m::Real)
Generate the spin state: $|j, m\rangle$
The eigenstate of the Spin-j
$\hat{S}_z$ operator with eigenvalue m
, where where j
is the spin quantum number and can be a non-negative integer or half-integer
See also jmat
.
sourceQuantumToolbox.spin_coherent
— Functionspin_coherent(j::Real, θ::Real, ϕ::Real)
Generate the coherent spin state (rotation of the $|j, j\rangle$ state), namely
\[|\theta, \phi \rangle = \hat{R}(\theta, \phi) |j, j\rangle\]
where the rotation operator is defined as
\[\hat{R}(\theta, \phi) = \exp \left( \frac{\theta}{2} (\hat{S}_- e^{i\phi} - \hat{S}_+ e^{-i\phi}) \right)\]
and $\hat{S}_\pm$ are plus and minus Spin-j
operators, respectively.
Arguments
j::Real
: The spin quantum number and can be a non-negative integer or half-integerθ::Real
: rotation angle from z-axisϕ::Real
: rotation angle from x-axis
See also jmat
and spin_state
.
Reference
sourceQuantumToolbox.bell_state
— Functionbell_state(x::Union{Int}, z::Union{Int})
Return the Bell state depending on the arguments (x, z)
:
(0, 0)
: $| \Phi^+ \rangle = ( |00\rangle + |11\rangle ) / \sqrt{2}$(0, 1)
: $| \Phi^- \rangle = ( |00\rangle - |11\rangle ) / \sqrt{2}$(1, 0)
: $| \Psi^+ \rangle = ( |01\rangle + |10\rangle ) / \sqrt{2}$(1, 1)
: $| \Psi^- \rangle = ( |01\rangle - |10\rangle ) / \sqrt{2}$
Here, x = 1
(z = 1
) means applying Pauli-$X$ ( Pauli-$Z$) unitary transformation on $| \Phi^+ \rangle$.
Example
julia> bell_state(0, 0)
Quantum Object: type=Ket dims=[2, 2] size=(4,)
4-element Vector{ComplexF64}:
0.7071067811865475 + 0.0im
@@ -196,7 +196,7 @@
0.0 + 0.0im
0.7071067811865475 + 0.0im
0.7071067811865475 + 0.0im
- 0.0 + 0.0im
Beware of type-stability! If you want to keep type stability, it is recommended to use bell_state(Val(x), Val(z))
instead of bell_state(x, z)
. See this link and the related Section for more details.
sourceQuantumToolbox.singlet_state
— Functionsinglet_state()
Return the two particle singlet state: $\frac{1}{\sqrt{2}} ( |01\rangle - |10\rangle )$
sourceQuantumToolbox.triplet_states
— Functiontriplet_states()
Return a list of the two particle triplet states:
- $|11\rangle$
- $( |01\rangle + |10\rangle ) / \sqrt{2}$
- $|00\rangle$
sourceQuantumToolbox.w_state
— Functionw_state(n::Union{Int,Val})
Returns the n
-qubit W-state:
\[\frac{1}{\sqrt{n}} \left( |100...0\rangle + |010...0\rangle + \cdots + |00...01\rangle \right)\]
Beware of type-stability! If you want to keep type stability, it is recommended to use w_state(Val(n))
instead of w_state(n)
. See this link and the related Section for more details.
sourceQuantumToolbox.ghz_state
— Functionghz_state(n::Union{Int,Val}; d::Int=2)
Returns the generalized n
-qudit Greenberger–Horne–Zeilinger (GHZ) state:
\[\frac{1}{\sqrt{d}} \sum_{i=0}^{d-1} | i \rangle \otimes \cdots \otimes | i \rangle\]
Here, d
specifies the dimension of each qudit. Default to d=2
(qubit).
Beware of type-stability! If you want to keep type stability, it is recommended to use ghz_state(Val(n))
instead of ghz_state(n)
. See this link and the related Section for more details.
sourceQuantumToolbox.rand_unitary
— Functionrand_unitary(dimensions, distribution=Val(:haar))
Returns a random unitary QuantumObject
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
The distribution
specifies which of the method used to obtain the unitary matrix:
:haar
: Haar random unitary matrix using the algorithm from reference 1:exp
: Uses $\exp(-i\hat{H})$, where $\hat{H}$ is a randomly generated Hermitian operator.
References
Beware of type-stability! If you want to keep type stability, it is recommended to use rand_unitary(dimensions, Val(distribution))
instead of rand_unitary(dimensions, distribution)
. Also, put dimensions
as Tuple
or SVector
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.sigmap
— Functionsigmap()
Pauli ladder operator $\hat{\sigma}_+ = (\hat{\sigma}_x + i \hat{\sigma}_y) / 2$.
See also jmat
.
sourceQuantumToolbox.sigmam
— Functionsigmam()
Pauli ladder operator $\hat{\sigma}_- = (\hat{\sigma}_x - i \hat{\sigma}_y) / 2$.
See also jmat
.
sourceQuantumToolbox.sigmax
— Functionsource QuantumToolbox.sigmay
— Functionsigmay()
Pauli operator $\hat{\sigma}_y = i \left( \hat{\sigma}_- - \hat{\sigma}_+ \right)$.
See also jmat
.
sourceQuantumToolbox.sigmaz
— Functionsource QuantumToolbox.jmat
— Functionjmat(j::Real, which::Union{Symbol,Val})
Generate higher-order Spin-j
operators, where j
is the spin quantum number and can be a non-negative integer or half-integer
The parameter which
specifies which of the following operator to return.
:x
: $\hat{S}_x$:y
: $\hat{S}_y$:z
: $\hat{S}_z$:+
: $\hat{S}_+$:-
: $\hat{S}_-$
Note that if the parameter which
is not specified, returns a set of Spin-j
operators: $(\hat{S}_x, \hat{S}_y, \hat{S}_z)$
Examples
julia> jmat(0.5, :x)
+ 0.0 + 0.0im
Beware of type-stability! If you want to keep type stability, it is recommended to use bell_state(Val(x), Val(z))
instead of bell_state(x, z)
. See this link and the related Section for more details.
sourceQuantumToolbox.singlet_state
— Functionsinglet_state()
Return the two particle singlet state: $\frac{1}{\sqrt{2}} ( |01\rangle - |10\rangle )$
sourceQuantumToolbox.triplet_states
— Functiontriplet_states()
Return a list of the two particle triplet states:
- $|11\rangle$
- $( |01\rangle + |10\rangle ) / \sqrt{2}$
- $|00\rangle$
sourceQuantumToolbox.w_state
— Functionw_state(n::Union{Int,Val})
Returns the n
-qubit W-state:
\[\frac{1}{\sqrt{n}} \left( |100...0\rangle + |010...0\rangle + \cdots + |00...01\rangle \right)\]
Beware of type-stability! If you want to keep type stability, it is recommended to use w_state(Val(n))
instead of w_state(n)
. See this link and the related Section for more details.
sourceQuantumToolbox.ghz_state
— Functionghz_state(n::Union{Int,Val}; d::Int=2)
Returns the generalized n
-qudit Greenberger–Horne–Zeilinger (GHZ) state:
\[\frac{1}{\sqrt{d}} \sum_{i=0}^{d-1} | i \rangle \otimes \cdots \otimes | i \rangle\]
Here, d
specifies the dimension of each qudit. Default to d=2
(qubit).
Beware of type-stability! If you want to keep type stability, it is recommended to use ghz_state(Val(n))
instead of ghz_state(n)
. See this link and the related Section for more details.
sourceQuantumToolbox.rand_unitary
— Functionrand_unitary(dimensions, distribution=Val(:haar))
Returns a random unitary QuantumObject
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
The distribution
specifies which of the method used to obtain the unitary matrix:
:haar
: Haar random unitary matrix using the algorithm from reference 1:exp
: Uses $\exp(-i\hat{H})$, where $\hat{H}$ is a randomly generated Hermitian operator.
References
Beware of type-stability! If you want to keep type stability, it is recommended to use rand_unitary(dimensions, Val(distribution))
instead of rand_unitary(dimensions, distribution)
. Also, put dimensions
as Tuple
or SVector
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.sigmap
— Functionsigmap()
Pauli ladder operator $\hat{\sigma}_+ = (\hat{\sigma}_x + i \hat{\sigma}_y) / 2$.
See also jmat
.
sourceQuantumToolbox.sigmam
— Functionsigmam()
Pauli ladder operator $\hat{\sigma}_- = (\hat{\sigma}_x - i \hat{\sigma}_y) / 2$.
See also jmat
.
sourceQuantumToolbox.sigmax
— Functionsource QuantumToolbox.sigmay
— Functionsigmay()
Pauli operator $\hat{\sigma}_y = i \left( \hat{\sigma}_- - \hat{\sigma}_+ \right)$.
See also jmat
.
sourceQuantumToolbox.sigmaz
— Functionsource QuantumToolbox.jmat
— Functionjmat(j::Real, which::Union{Symbol,Val})
Generate higher-order Spin-j
operators, where j
is the spin quantum number and can be a non-negative integer or half-integer
The parameter which
specifies which of the following operator to return.
:x
: $\hat{S}_x$:y
: $\hat{S}_y$:z
: $\hat{S}_z$:+
: $\hat{S}_+$:-
: $\hat{S}_-$
Note that if the parameter which
is not specified, returns a set of Spin-j
operators: $(\hat{S}_x, \hat{S}_y, \hat{S}_z)$
Examples
julia> jmat(0.5, :x)
Quantum Object: type=Operator dims=[2] size=(2, 2) ishermitian=true
2×2 SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
⋅ 0.5+0.0im
@@ -214,7 +214,7 @@
1.5+0.0im ⋅ ⋅ ⋅
⋅ 0.5+0.0im ⋅ ⋅
⋅ ⋅ -0.5+0.0im ⋅
- ⋅ ⋅ ⋅ -1.5+0.0im
Beware of type-stability! If you want to keep type stability, it is recommended to use jmat(j, Val(which))
instead of jmat(j, which)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.spin_Jx
— Functionspin_Jx(j::Real)
$\hat{S}_x$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_Jy
— Functionspin_Jy(j::Real)
$\hat{S}_y$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_Jz
— Functionspin_Jz(j::Real)
$\hat{S}_z$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_Jm
— Functionspin_Jm(j::Real)
$\hat{S}_-$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_Jp
— Functionspin_Jp(j::Real)
$\hat{S}_+$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_J_set
— Functionspin_J_set(j::Real)
A set of Spin-j
operators $(\hat{S}_x, \hat{S}_y, \hat{S}_z)$, where j
is the spin quantum number and can be a non-negative integer or half-integer.
Note that this functions is same as jmat(j)
. See also jmat
.
sourceQuantumToolbox.destroy
— Functiondestroy(N::Int)
Bosonic annihilation operator with Hilbert space cutoff N
.
This operator acts on a fock state as $\hat{a} \ket{n} = \sqrt{n} \ket{n-1}$.
Examples
julia> a = destroy(20)
+ ⋅ ⋅ ⋅ -1.5+0.0im
Beware of type-stability! If you want to keep type stability, it is recommended to use jmat(j, Val(which))
instead of jmat(j, which)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.spin_Jx
— Functionspin_Jx(j::Real)
$\hat{S}_x$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_Jy
— Functionspin_Jy(j::Real)
$\hat{S}_y$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_Jz
— Functionspin_Jz(j::Real)
$\hat{S}_z$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_Jm
— Functionspin_Jm(j::Real)
$\hat{S}_-$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_Jp
— Functionspin_Jp(j::Real)
$\hat{S}_+$ operator for Spin-j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
sourceQuantumToolbox.spin_J_set
— Functionspin_J_set(j::Real)
A set of Spin-j
operators $(\hat{S}_x, \hat{S}_y, \hat{S}_z)$, where j
is the spin quantum number and can be a non-negative integer or half-integer.
Note that this functions is same as jmat(j)
. See also jmat
.
sourceQuantumToolbox.destroy
— Functiondestroy(N::Int)
Bosonic annihilation operator with Hilbert space cutoff N
.
This operator acts on a fock state as $\hat{a} \ket{n} = \sqrt{n} \ket{n-1}$.
Examples
julia> a = destroy(20)
Quantum Object: type=Operator dims=[20] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 19 stored entries:
⎡⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀⎤
@@ -224,7 +224,7 @@
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢⎦
julia> fock(20, 3)' * a * fock(20, 4)
-2.0 + 0.0im
sourceQuantumToolbox.create
— Functioncreate(N::Int)
Bosonic creation operator with Hilbert space cutoff N
.
This operator acts on a fock state as $\hat{a}^\dagger \ket{n} = \sqrt{n+1} \ket{n+1}$.
Examples
julia> a_d = create(20)
+2.0 + 0.0im
sourceQuantumToolbox.create
— Functioncreate(N::Int)
Bosonic creation operator with Hilbert space cutoff N
.
This operator acts on a fock state as $\hat{a}^\dagger \ket{n} = \sqrt{n+1} \ket{n+1}$.
Examples
julia> a_d = create(20)
Quantum Object: type=Operator dims=[20] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 19 stored entries:
⎡⠢⡀⠀⠀⠀⠀⠀⠀⠀⠀⎤
@@ -234,15 +234,15 @@
⎣⠀⠀⠀⠀⠀⠀⠀⠈⠢⡀⎦
julia> fock(20, 4)' * a_d * fock(20, 3)
-2.0 + 0.0im
sourceQuantumToolbox.displace
— Functiondisplace(N::Int, α::Number)
Generate a displacement operator:
\[\hat{D}(\alpha)=\exp\left( \alpha \hat{a}^\dagger - \alpha^* \hat{a} \right),\]
where $\hat{a}$ is the bosonic annihilation operator, and $\alpha$ is the amount of displacement in optical phase space.
sourceQuantumToolbox.squeeze
— Functionsqueeze(N::Int, z::Number)
Generate a single-mode squeeze operator:
\[\hat{S}(z)=\exp\left( \frac{1}{2} (z^* \hat{a}^2 - z(\hat{a}^\dagger)^2) \right),\]
where $\hat{a}$ is the bosonic annihilation operator.
sourceQuantumToolbox.num
— Functionnum(N::Int)
Bosonic number operator with Hilbert space cutoff N
.
This operator is defined as $\hat{N}=\hat{a}^\dagger \hat{a}$, where $\hat{a}$ is the bosonic annihilation operator.
sourceQuantumToolbox.position
— Functionposition(N::Int)
Position operator with Hilbert space cutoff N
.
This operator is defined as $\hat{x}=\frac{1}{\sqrt{2}} (\hat{a}^\dagger + \hat{a})$, where $\hat{a}$ is the bosonic annihilation operator.
sourceQuantumToolbox.momentum
— Functionmomentum(N::Int)
Momentum operator with Hilbert space cutoff N
.
This operator is defined as $\hat{p}= \frac{i}{\sqrt{2}} (\hat{a}^\dagger - \hat{a})$, where $\hat{a}$ is the bosonic annihilation operator.
sourceQuantumToolbox.phase
— Functionphase(N::Int, ϕ0::Real=0)
Single-mode Pegg-Barnett phase operator with Hilbert space cutoff $N$ and the reference phase $\phi_0$.
This operator is defined as
\[\hat{\phi} = \sum_{m=0}^{N-1} \phi_m |\phi_m\rangle \langle\phi_m|,\]
where
\[\phi_m = \phi_0 + \frac{2m\pi}{N},\]
and
\[|\phi_m\rangle = \frac{1}{\sqrt{N}} \sum_{n=0}^{N-1} \exp(i n \phi_m) |n\rangle.\]
Reference
sourceQuantumToolbox.fdestroy
— Functionfdestroy(N::Union{Int,Val}, j::Int)
Construct a fermionic destruction operator acting on the j
-th site, where the fock space has totally N
-sites:
Here, we use the Jordan-Wigner transformation, namely
\[\hat{d}_j = \hat{\sigma}_z^{\otimes j-1} \otimes \hat{\sigma}_{+} \otimes \hat{\mathbb{1}}^{\otimes N-j}\]
The site index j
should satisfy: 1 ≤ j ≤ N
.
Note that we put $\hat{\sigma}_{+} = \begin{pmatrix} 0 & 1 \\ 0 & 0 \end{pmatrix}$ here because we consider $|0\rangle = \begin{pmatrix} 1 \\ 0 \end{pmatrix}$ to be ground (vacant) state, and $|1\rangle = \begin{pmatrix} 0 \\ 1 \end{pmatrix}$ to be excited (occupied) state.
Beware of type-stability! If you want to keep type stability, it is recommended to use fdestroy(Val(N), j)
instead of fdestroy(N, j)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.fcreate
— Functionfcreate(N::Union{Int,Val}, j::Int)
Construct a fermionic creation operator acting on the j
-th site, where the fock space has totally N
-sites:
Here, we use the Jordan-Wigner transformation, namely
\[\hat{d}^\dagger_j = \hat{\sigma}_z^{\otimes j-1} \otimes \hat{\sigma}_{-} \otimes \hat{\mathbb{1}}^{\otimes N-j}\]
The site index j
should satisfy: 1 ≤ j ≤ N
.
Note that we put $\hat{\sigma}_{-} = \begin{pmatrix} 0 & 0 \\ 1 & 0 \end{pmatrix}$ here because we consider $|0\rangle = \begin{pmatrix} 1 \\ 0 \end{pmatrix}$ to be ground (vacant) state, and $|1\rangle = \begin{pmatrix} 0 \\ 1 \end{pmatrix}$ to be excited (occupied) state.
Beware of type-stability! If you want to keep type stability, it is recommended to use fcreate(Val(N), j)
instead of fcreate(N, j)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.tunneling
— Functiontunneling(N::Int, m::Int=1; sparse::Union{Bool,Val{<:Bool}}=Val(false))
Generate a tunneling operator defined as:
\[\sum_{n=0}^{N-m} | n \rangle\langle n+m | + | n+m \rangle\langle n |,\]
where $N$ is the number of basis states in the Hilbert space, and $m$ is the number of excitations in tunneling event.
If sparse=true
, the operator is returned as a sparse matrix, otherwise a dense matrix is returned.
Beware of type-stability! If you want to keep type stability, it is recommended to use tunneling(N, m, Val(sparse))
instead of tunneling(N, m, sparse)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.qft
— Functionqft(dimensions)
Generates a discrete Fourier transform matrix $\hat{F}_N$ for Quantum Fourier Transform (QFT) with given argument dimensions
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
$N$ represents the total dimension, and therefore the matrix is defined as
\[\hat{F}_N = \frac{1}{\sqrt{N}}\begin{bmatrix}
+2.0 + 0.0im
sourceQuantumToolbox.displace
— Functiondisplace(N::Int, α::Number)
Generate a displacement operator:
\[\hat{D}(\alpha)=\exp\left( \alpha \hat{a}^\dagger - \alpha^* \hat{a} \right),\]
where $\hat{a}$ is the bosonic annihilation operator, and $\alpha$ is the amount of displacement in optical phase space.
sourceQuantumToolbox.squeeze
— Functionsqueeze(N::Int, z::Number)
Generate a single-mode squeeze operator:
\[\hat{S}(z)=\exp\left( \frac{1}{2} (z^* \hat{a}^2 - z(\hat{a}^\dagger)^2) \right),\]
where $\hat{a}$ is the bosonic annihilation operator.
sourceQuantumToolbox.num
— Functionnum(N::Int)
Bosonic number operator with Hilbert space cutoff N
.
This operator is defined as $\hat{N}=\hat{a}^\dagger \hat{a}$, where $\hat{a}$ is the bosonic annihilation operator.
sourceQuantumToolbox.position
— Functionposition(N::Int)
Position operator with Hilbert space cutoff N
.
This operator is defined as $\hat{x}=\frac{1}{\sqrt{2}} (\hat{a}^\dagger + \hat{a})$, where $\hat{a}$ is the bosonic annihilation operator.
sourceQuantumToolbox.momentum
— Functionmomentum(N::Int)
Momentum operator with Hilbert space cutoff N
.
This operator is defined as $\hat{p}= \frac{i}{\sqrt{2}} (\hat{a}^\dagger - \hat{a})$, where $\hat{a}$ is the bosonic annihilation operator.
sourceQuantumToolbox.phase
— Functionphase(N::Int, ϕ0::Real=0)
Single-mode Pegg-Barnett phase operator with Hilbert space cutoff $N$ and the reference phase $\phi_0$.
This operator is defined as
\[\hat{\phi} = \sum_{m=0}^{N-1} \phi_m |\phi_m\rangle \langle\phi_m|,\]
where
\[\phi_m = \phi_0 + \frac{2m\pi}{N},\]
and
\[|\phi_m\rangle = \frac{1}{\sqrt{N}} \sum_{n=0}^{N-1} \exp(i n \phi_m) |n\rangle.\]
Reference
sourceQuantumToolbox.fdestroy
— Functionfdestroy(N::Union{Int,Val}, j::Int)
Construct a fermionic destruction operator acting on the j
-th site, where the fock space has totally N
-sites:
Here, we use the Jordan-Wigner transformation, namely
\[\hat{d}_j = \hat{\sigma}_z^{\otimes j-1} \otimes \hat{\sigma}_{+} \otimes \hat{\mathbb{1}}^{\otimes N-j}\]
The site index j
should satisfy: 1 ≤ j ≤ N
.
Note that we put $\hat{\sigma}_{+} = \begin{pmatrix} 0 & 1 \\ 0 & 0 \end{pmatrix}$ here because we consider $|0\rangle = \begin{pmatrix} 1 \\ 0 \end{pmatrix}$ to be ground (vacant) state, and $|1\rangle = \begin{pmatrix} 0 \\ 1 \end{pmatrix}$ to be excited (occupied) state.
Beware of type-stability! If you want to keep type stability, it is recommended to use fdestroy(Val(N), j)
instead of fdestroy(N, j)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.fcreate
— Functionfcreate(N::Union{Int,Val}, j::Int)
Construct a fermionic creation operator acting on the j
-th site, where the fock space has totally N
-sites:
Here, we use the Jordan-Wigner transformation, namely
\[\hat{d}^\dagger_j = \hat{\sigma}_z^{\otimes j-1} \otimes \hat{\sigma}_{-} \otimes \hat{\mathbb{1}}^{\otimes N-j}\]
The site index j
should satisfy: 1 ≤ j ≤ N
.
Note that we put $\hat{\sigma}_{-} = \begin{pmatrix} 0 & 0 \\ 1 & 0 \end{pmatrix}$ here because we consider $|0\rangle = \begin{pmatrix} 1 \\ 0 \end{pmatrix}$ to be ground (vacant) state, and $|1\rangle = \begin{pmatrix} 0 \\ 1 \end{pmatrix}$ to be excited (occupied) state.
Beware of type-stability! If you want to keep type stability, it is recommended to use fcreate(Val(N), j)
instead of fcreate(N, j)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.tunneling
— Functiontunneling(N::Int, m::Int=1; sparse::Union{Bool,Val{<:Bool}}=Val(false))
Generate a tunneling operator defined as:
\[\sum_{n=0}^{N-m} | n \rangle\langle n+m | + | n+m \rangle\langle n |,\]
where $N$ is the number of basis states in the Hilbert space, and $m$ is the number of excitations in tunneling event.
If sparse=true
, the operator is returned as a sparse matrix, otherwise a dense matrix is returned.
Beware of type-stability! If you want to keep type stability, it is recommended to use tunneling(N, m, Val(sparse))
instead of tunneling(N, m, sparse)
. See this link and the related Section about type stability for more details.
sourceQuantumToolbox.qft
— Functionqft(dimensions)
Generates a discrete Fourier transform matrix $\hat{F}_N$ for Quantum Fourier Transform (QFT) with given argument dimensions
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
$N$ represents the total dimension, and therefore the matrix is defined as
\[\hat{F}_N = \frac{1}{\sqrt{N}}\begin{bmatrix}
1 & 1 & 1 & 1 & \cdots & 1\\
1 & \omega & \omega^2 & \omega^3 & \cdots & \omega^{N-1}\\
1 & \omega^2 & \omega^4 & \omega^6 & \cdots & \omega^{2(N-1)}\\
1 & \omega^3 & \omega^6 & \omega^9 & \cdots & \omega^{3(N-1)}\\
\vdots & \vdots & \vdots & \vdots & \ddots & \vdots\\
1 & \omega^{N-1} & \omega^{2(N-1)} & \omega^{3(N-1)} & \cdots & \omega^{(N-1)(N-1)}
-\end{bmatrix},\]
where $\omega = \exp(\frac{2 \pi i}{N})$.
Beware of type-stability! It is highly recommended to use qft(dimensions)
with dimensions
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.eye
— Functioneye(N::Int; type=Operator, dims=nothing)
Identity operator $\hat{\mathbb{1}}$ with size N
.
It is also possible to specify the list of Hilbert dimensions dims
if different subsystems are present.
Note that type
can only be either Operator
or SuperOperator
sourceQuantumToolbox.projection
— Functionprojection(N::Int, i::Int, j::Int)
Generates the projection operator $\hat{O} = \dyad{i}{j}$ with Hilbert space dimension N
.
sourceQuantumToolbox.commutator
— Functioncommutator(A::QuantumObject, B::QuantumObject; anti::Bool=false)
Return the commutator (or anti
-commutator) of the two QuantumObject
:
- commutator (
anti=false
): $\hat{A}\hat{B}-\hat{B}\hat{A}$ - anticommutator (
anti=true
): $\hat{A}\hat{B}+\hat{B}\hat{A}$
Note that A
and B
must be Operator
sourceQuantumToolbox.spre
— Functionspre(A::QuantumObject, Id_cache=I(size(A,1)))
Returns the SuperOperator
form of A
acting on the left of the density matrix operator: $\mathcal{O} \left(\hat{A}\right) \left[ \hat{\rho} \right] = \hat{A} \hat{\rho}$.
Since the density matrix is vectorized in OperatorKet
form: $|\hat{\rho}\rangle\rangle$, this SuperOperator
is always a matrix $\hat{\mathbb{1}} \otimes \hat{A}$, namely
\[\mathcal{O} \left(\hat{A}\right) \left[ \hat{\rho} \right] = \hat{\mathbb{1}} \otimes \hat{A} ~ |\hat{\rho}\rangle\rangle\]
(see the section in documentation: Superoperators and Vectorized Operators for more details)
The optional argument Id_cache
can be used to pass a precomputed identity matrix. This can be useful when the same function is applied multiple times with a known Hilbert space dimension.
sourceQuantumToolbox.spost
— Functionspost(B::QuantumObject, Id_cache=I(size(B,1)))
Returns the SuperOperator
form of B
acting on the right of the density matrix operator: $\mathcal{O} \left(\hat{B}\right) \left[ \hat{\rho} \right] = \hat{\rho} \hat{B}$.
Since the density matrix is vectorized in OperatorKet
form: $|\hat{\rho}\rangle\rangle$, this SuperOperator
is always a matrix $\hat{B}^T \otimes \hat{\mathbb{1}}$, namely
\[\mathcal{O} \left(\hat{B}\right) \left[ \hat{\rho} \right] = \hat{B}^T \otimes \hat{\mathbb{1}} ~ |\hat{\rho}\rangle\rangle\]
(see the section in documentation: Superoperators and Vectorized Operators for more details)
The optional argument Id_cache
can be used to pass a precomputed identity matrix. This can be useful when the same function is applied multiple times with a known Hilbert space dimension.
sourceQuantumToolbox.sprepost
— Functionsprepost(A::QuantumObject, B::QuantumObject)
Returns the SuperOperator
form of A
and B
acting on the left and right of the density matrix operator, respectively: $\mathcal{O} \left( \hat{A}, \hat{B} \right) \left[ \hat{\rho} \right] = \hat{A} \hat{\rho} \hat{B}$.
Since the density matrix is vectorized in OperatorKet
form: $|\hat{\rho}\rangle\rangle$, this SuperOperator
is always a matrix $\hat{B}^T \otimes \hat{A}$, namely
\[\mathcal{O} \left(\hat{A}, \hat{B}\right) \left[ \hat{\rho} \right] = \hat{B}^T \otimes \hat{A} ~ |\hat{\rho}\rangle\rangle = \textrm{spre}(\hat{A}) * \textrm{spost}(\hat{B}) ~ |\hat{\rho}\rangle\rangle\]
(see the section in documentation: Superoperators and Vectorized Operators for more details)
sourceQuantumToolbox.lindblad_dissipator
— Functionlindblad_dissipator(O::QuantumObject, Id_cache=I(size(O,1))
Returns the Lindblad SuperOperator
defined as
\[\mathcal{D} \left( \hat{O} \right) \left[ \hat{\rho} \right] = \frac{1}{2} \left( 2 \hat{O} \hat{\rho} \hat{O}^\dagger -
-\hat{O}^\dagger \hat{O} \hat{\rho} - \hat{\rho} \hat{O}^\dagger \hat{O} \right)\]
The optional argument Id_cache
can be used to pass a precomputed identity matrix. This can be useful when the same function is applied multiple times with a known Hilbert space dimension.
sourceSynonyms of functions for Qobj
QuantumToolbox.Qobj
— FunctionQobj(A::AbstractArray; type::QuantumObjectType, dims::Vector{Int})
Generate QuantumObject
Note that this functions is same as QuantumObject(A; type=type, dims=dims)
sourceQuantumToolbox.shape
— Functionshape(A::QuantumObject)
Returns a tuple containing each dimensions of the array in the QuantumObject
.
Note that this function is same as size(A)
sourceQuantumToolbox.isherm
— Functionisherm(A::QuantumObject)
Test whether the QuantumObject
is Hermitian.
Note that this functions is same as ishermitian(A)
sourceQuantumToolbox.trans
— Functiontrans(A::QuantumObject)
Lazy matrix transpose of the QuantumObject
.
Note that this function is same as transpose(A)
sourceQuantumToolbox.dag
— Functiondag(A::QuantumObject)
Lazy adjoint (conjugate transposition) of the QuantumObject
Note that this function is same as adjoint(A)
sourceQuantumToolbox.matrix_element
— Functionmatrix_element(i::QuantumObject, A::QuantumObject j::QuantumObject)
Compute the generalized dot product dot(i, A*j)
between three QuantumObject
: $\langle i | \hat{A} | j \rangle$
Note that this function is same as dot(i, A, j)
Supports the following inputs:
A
is in the type of Operator
, with i
and j
are both Ket
.A
is in the type of SuperOperator
, with i
and j
are both OperatorKet
sourceQuantumToolbox.unit
— Functionunit(A::QuantumObject, p::Real)
Return normalized QuantumObject
so that its p
-norm equals to unity, i.e. norm(A, p) == 1
.
Support for the following types of QuantumObject
:
Note that this function is same as normalize(A, p)
Also, see norm
about its definition for different types of QuantumObject
.
sourceQuantumToolbox.sqrtm
— Functionsqrtm(A::QuantumObject)
Matrix square root of Operator
type of QuantumObject
Note that for other types of QuantumObject
use sprt(A)
instead.
sourceQuantumToolbox.logm
— Functionlogm(A::QuantumObject)
Matrix logarithm of QuantumObject
Note that this function is same as log(A)
and only supports for Operator
and SuperOperator
sourceQuantumToolbox.expm
— Functionexpm(A::QuantumObject)
Matrix exponential of QuantumObject
Note that this function is same as exp(A)
and only supports for Operator
and SuperOperator
sourceQuantumToolbox.sinm
— Functionsinm(A::QuantumObject)
Matrix sine of QuantumObject
, defined as
$\sin \left( \hat{A} \right) = \frac{e^{i \hat{A}} - e^{-i \hat{A}}}{2 i}$
Note that this function is same as sin(A)
and only supports for Operator
and SuperOperator
sourceQuantumToolbox.cosm
— Functioncosm(A::QuantumObject)
Matrix cosine of QuantumObject
, defined as
$\cos \left( \hat{A} \right) = \frac{e^{i \hat{A}} + e^{-i \hat{A}}}{2}$
Note that this function is same as cos(A)
and only supports for Operator
and SuperOperator
sourceQuantumToolbox.tensor
— Functiontensor(A::QuantumObject, B::QuantumObject, ...)
Returns the Kronecker product $\hat{A} \otimes \hat{B} \otimes \cdots$.
Note that this function is same as kron(A, B, ...)
Examples
julia> x = sigmax()
+\end{bmatrix},\]where $\omega = \exp(\frac{2 \pi i}{N})$.
Beware of type-stability! It is highly recommended to use qft(dimensions)
with dimensions
as Tuple
or SVector
to keep type stability. See the related Section about type stability for more details.
sourceQuantumToolbox.eye
— Functioneye(N::Int; type=Operator, dims=nothing)
Identity operator $\hat{\mathbb{1}}$ with size N
.
It is also possible to specify the list of Hilbert dimensions dims
if different subsystems are present.
Note that type
can only be either Operator
or SuperOperator
sourceQuantumToolbox.projection
— Functionprojection(N::Int, i::Int, j::Int)
Generates the projection operator $\hat{O} = \dyad{i}{j}$ with Hilbert space dimension N
.
sourceQuantumToolbox.commutator
— Functioncommutator(A::QuantumObject, B::QuantumObject; anti::Bool=false)
Return the commutator (or anti
-commutator) of the two QuantumObject
:
- commutator (
anti=false
): $\hat{A}\hat{B}-\hat{B}\hat{A}$ - anticommutator (
anti=true
): $\hat{A}\hat{B}+\hat{B}\hat{A}$
Note that A
and B
must be Operator
sourceQuantumToolbox.spre
— Functionspre(A::QuantumObject, Id_cache=I(size(A,1)))
Returns the SuperOperator
form of A
acting on the left of the density matrix operator: $\mathcal{O} \left(\hat{A}\right) \left[ \hat{\rho} \right] = \hat{A} \hat{\rho}$.
Since the density matrix is vectorized in OperatorKet
form: $|\hat{\rho}\rangle\rangle$, this SuperOperator
is always a matrix $\hat{\mathbb{1}} \otimes \hat{A}$, namely
\[\mathcal{O} \left(\hat{A}\right) \left[ \hat{\rho} \right] = \hat{\mathbb{1}} \otimes \hat{A} ~ |\hat{\rho}\rangle\rangle\]
(see the section in documentation: Superoperators and Vectorized Operators for more details)
The optional argument Id_cache
can be used to pass a precomputed identity matrix. This can be useful when the same function is applied multiple times with a known Hilbert space dimension.
sourceQuantumToolbox.spost
— Functionspost(B::QuantumObject, Id_cache=I(size(B,1)))
Returns the SuperOperator
form of B
acting on the right of the density matrix operator: $\mathcal{O} \left(\hat{B}\right) \left[ \hat{\rho} \right] = \hat{\rho} \hat{B}$.
Since the density matrix is vectorized in OperatorKet
form: $|\hat{\rho}\rangle\rangle$, this SuperOperator
is always a matrix $\hat{B}^T \otimes \hat{\mathbb{1}}$, namely
\[\mathcal{O} \left(\hat{B}\right) \left[ \hat{\rho} \right] = \hat{B}^T \otimes \hat{\mathbb{1}} ~ |\hat{\rho}\rangle\rangle\]
(see the section in documentation: Superoperators and Vectorized Operators for more details)
The optional argument Id_cache
can be used to pass a precomputed identity matrix. This can be useful when the same function is applied multiple times with a known Hilbert space dimension.
sourceQuantumToolbox.sprepost
— Functionsprepost(A::QuantumObject, B::QuantumObject)
Returns the SuperOperator
form of A
and B
acting on the left and right of the density matrix operator, respectively: $\mathcal{O} \left( \hat{A}, \hat{B} \right) \left[ \hat{\rho} \right] = \hat{A} \hat{\rho} \hat{B}$.
Since the density matrix is vectorized in OperatorKet
form: $|\hat{\rho}\rangle\rangle$, this SuperOperator
is always a matrix $\hat{B}^T \otimes \hat{A}$, namely
\[\mathcal{O} \left(\hat{A}, \hat{B}\right) \left[ \hat{\rho} \right] = \hat{B}^T \otimes \hat{A} ~ |\hat{\rho}\rangle\rangle = \textrm{spre}(\hat{A}) * \textrm{spost}(\hat{B}) ~ |\hat{\rho}\rangle\rangle\]
(see the section in documentation: Superoperators and Vectorized Operators for more details)
sourceQuantumToolbox.lindblad_dissipator
— Functionlindblad_dissipator(O::QuantumObject, Id_cache=I(size(O,1))
Returns the Lindblad SuperOperator
defined as
\[\mathcal{D} \left( \hat{O} \right) \left[ \hat{\rho} \right] = \frac{1}{2} \left( 2 \hat{O} \hat{\rho} \hat{O}^\dagger -
+\hat{O}^\dagger \hat{O} \hat{\rho} - \hat{\rho} \hat{O}^\dagger \hat{O} \right)\]
The optional argument Id_cache
can be used to pass a precomputed identity matrix. This can be useful when the same function is applied multiple times with a known Hilbert space dimension.
sourceSynonyms of functions for Qobj
QuantumToolbox.Qobj
— FunctionQobj(A::AbstractArray; type::QuantumObjectType, dims::Vector{Int})
Generate QuantumObject
Note that this functions is same as QuantumObject(A; type=type, dims=dims)
sourceQuantumToolbox.shape
— Functionshape(A::QuantumObject)
Returns a tuple containing each dimensions of the array in the QuantumObject
.
Note that this function is same as size(A)
sourceQuantumToolbox.isherm
— Functionisherm(A::QuantumObject)
Test whether the QuantumObject
is Hermitian.
Note that this functions is same as ishermitian(A)
sourceQuantumToolbox.trans
— Functiontrans(A::QuantumObject)
Lazy matrix transpose of the QuantumObject
.
Note that this function is same as transpose(A)
sourceQuantumToolbox.dag
— Functiondag(A::QuantumObject)
Lazy adjoint (conjugate transposition) of the QuantumObject
Note that this function is same as adjoint(A)
sourceQuantumToolbox.matrix_element
— Functionmatrix_element(i::QuantumObject, A::QuantumObject j::QuantumObject)
Compute the generalized dot product dot(i, A*j)
between three QuantumObject
: $\langle i | \hat{A} | j \rangle$
Note that this function is same as dot(i, A, j)
Supports the following inputs:
A
is in the type of Operator
, with i
and j
are both Ket
.A
is in the type of SuperOperator
, with i
and j
are both OperatorKet
sourceQuantumToolbox.unit
— Functionunit(A::QuantumObject, p::Real)
Return normalized QuantumObject
so that its p
-norm equals to unity, i.e. norm(A, p) == 1
.
Support for the following types of QuantumObject
:
Note that this function is same as normalize(A, p)
Also, see norm
about its definition for different types of QuantumObject
.
sourceQuantumToolbox.sqrtm
— Functionsqrtm(A::QuantumObject)
Matrix square root of Operator
type of QuantumObject
Note that for other types of QuantumObject
use sprt(A)
instead.
sourceQuantumToolbox.logm
— Functionlogm(A::QuantumObject)
Matrix logarithm of QuantumObject
Note that this function is same as log(A)
and only supports for Operator
and SuperOperator
sourceQuantumToolbox.expm
— Functionexpm(A::QuantumObject)
Matrix exponential of QuantumObject
Note that this function is same as exp(A)
and only supports for Operator
and SuperOperator
sourceQuantumToolbox.sinm
— Functionsinm(A::QuantumObject)
Matrix sine of QuantumObject
, defined as
$\sin \left( \hat{A} \right) = \frac{e^{i \hat{A}} - e^{-i \hat{A}}}{2 i}$
Note that this function is same as sin(A)
and only supports for Operator
and SuperOperator
sourceQuantumToolbox.cosm
— Functioncosm(A::QuantumObject)
Matrix cosine of QuantumObject
, defined as
$\cos \left( \hat{A} \right) = \frac{e^{i \hat{A}} + e^{-i \hat{A}}}{2}$
Note that this function is same as cos(A)
and only supports for Operator
and SuperOperator
sourceQuantumToolbox.tensor
— Functiontensor(A::QuantumObject, B::QuantumObject, ...)
Returns the Kronecker product $\hat{A} \otimes \hat{B} \otimes \cdots$.
Note that this function is same as kron(A, B, ...)
Examples
julia> x = sigmax()
Quantum Object: type=Operator dims=[2] size=(2, 2) ishermitian=true
2×2 SparseMatrixCSC{ComplexF64, Int64} with 2 stored entries:
⋅ 1.0+0.0im
@@ -260,7 +260,7 @@
⋅ ⋅ ⋅ ⋅ ⋅ ⋅
⋅ ⋅ 1.0+0.0im … ⋅ ⋅ ⋅
⋅ 1.0+0.0im ⋅ ⋅ ⋅ ⋅
- 1.0+0.0im ⋅ ⋅ ⋅ ⋅ ⋅
sourceQuantumToolbox.:⊗
— Function⊗(A::QuantumObject, B::QuantumObject)
Returns the Kronecker product $\hat{A} \otimes \hat{B}$.
Note that this function is same as kron(A, B)
Examples
julia> a = destroy(20)
+ 1.0+0.0im ⋅ ⋅ ⋅ ⋅ ⋅
sourceQuantumToolbox.:⊗
— Function⊗(A::QuantumObject, B::QuantumObject)
Returns the Kronecker product $\hat{A} \otimes \hat{B}$.
Note that this function is same as kron(A, B)
Examples
julia> a = destroy(20)
Quantum Object: type=Operator dims=[20] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 19 stored entries:
⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀
@@ -292,7 +292,7 @@
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡀⠀
-⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⠦
sourceQuantumToolbox.qeye
— Functionqeye(N::Int; type=Operator, dims=nothing)
Identity operator $\hat{\mathbb{1}}$ with size N
.
It is also possible to specify the list of Hilbert dimensions dims
if different subsystems are present.
Note that this function is same as eye(N, type=type, dims=dims)
, and type
can only be either Operator
or SuperOperator
sourceTime evolution
QuantumToolbox.TimeEvolutionSol
— Typestruct TimeEvolutionSol
A structure storing the results and some information from solving time evolution.
Fields (Attributes)
times::AbstractVector
: The time list of the evolution.states::Vector{QuantumObject}
: The list of result states.expect::Matrix
: The expectation values corresponding to each time point in times
.retcode
: The return code from the solver.alg
: The algorithm which is used during the solving process.abstol::Real
: The absolute tolerance which is used during the solving process.reltol::Real
: The relative tolerance which is used during the solving process.
sourceQuantumToolbox.TimeEvolutionMCSol
— Typestruct TimeEvolutionMCSol
A structure storing the results and some information from solving quantum trajectories of the Monte Carlo wave function time evolution.
Fields (Attributes)
ntraj::Int
: Number of trajectoriestimes::AbstractVector
: The time list of the evolution.states::Vector{Vector{QuantumObject}}
: The list of result states in each trajectory.expect::Matrix
: The expectation values (averaging all trajectories) corresponding to each time point in times
.expect_all::Array
: The expectation values corresponding to each trajectory and each time point in times
jump_times::Vector{Vector{Real}}
: The time records of every quantum jump occurred in each trajectory.jump_which::Vector{Vector{Int}}
: The indices of the jump operators in c_ops
that describe the corresponding quantum jumps occurred in each trajectory.converged::Bool
: Whether the solution is converged or not.alg
: The algorithm which is used during the solving process.abstol::Real
: The absolute tolerance which is used during the solving process.reltol::Real
: The relative tolerance which is used during the solving process.
sourceQuantumToolbox.TimeEvolutionSSESol
— Typestruct TimeEvolutionSSESol
A structure storing the results and some information from solving trajectories of the Stochastic Shrodinger equation time evolution.
Fields (Attributes)
ntraj::Int
: Number of trajectoriestimes::AbstractVector
: The time list of the evolution.states::Vector{Vector{QuantumObject}}
: The list of result states in each trajectory.expect::Matrix
: The expectation values (averaging all trajectories) corresponding to each time point in times
.expect_all::Array
: The expectation values corresponding to each trajectory and each time point in times
converged::Bool
: Whether the solution is converged or not.alg
: The algorithm which is used during the solving process.abstol::Real
: The absolute tolerance which is used during the solving process.reltol::Real
: The relative tolerance which is used during the solving process.
sourceQuantumToolbox.sesolveProblem
— FunctionsesolveProblem(H::QuantumObject,
+⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⠦
sourceQuantumToolbox.qeye
— Functionqeye(N::Int; type=Operator, dims=nothing)
Identity operator $\hat{\mathbb{1}}$ with size N
.
It is also possible to specify the list of Hilbert dimensions dims
if different subsystems are present.
Note that this function is same as eye(N, type=type, dims=dims)
, and type
can only be either Operator
or SuperOperator
sourceTime evolution
QuantumToolbox.TimeEvolutionSol
— Typestruct TimeEvolutionSol
A structure storing the results and some information from solving time evolution.
Fields (Attributes)
times::AbstractVector
: The time list of the evolution.states::Vector{QuantumObject}
: The list of result states.expect::Matrix
: The expectation values corresponding to each time point in times
.retcode
: The return code from the solver.alg
: The algorithm which is used during the solving process.abstol::Real
: The absolute tolerance which is used during the solving process.reltol::Real
: The relative tolerance which is used during the solving process.
sourceQuantumToolbox.TimeEvolutionMCSol
— Typestruct TimeEvolutionMCSol
A structure storing the results and some information from solving quantum trajectories of the Monte Carlo wave function time evolution.
Fields (Attributes)
ntraj::Int
: Number of trajectoriestimes::AbstractVector
: The time list of the evolution.states::Vector{Vector{QuantumObject}}
: The list of result states in each trajectory.expect::Matrix
: The expectation values (averaging all trajectories) corresponding to each time point in times
.expect_all::Array
: The expectation values corresponding to each trajectory and each time point in times
jump_times::Vector{Vector{Real}}
: The time records of every quantum jump occurred in each trajectory.jump_which::Vector{Vector{Int}}
: The indices of the jump operators in c_ops
that describe the corresponding quantum jumps occurred in each trajectory.converged::Bool
: Whether the solution is converged or not.alg
: The algorithm which is used during the solving process.abstol::Real
: The absolute tolerance which is used during the solving process.reltol::Real
: The relative tolerance which is used during the solving process.
sourceQuantumToolbox.TimeEvolutionSSESol
— Typestruct TimeEvolutionSSESol
A structure storing the results and some information from solving trajectories of the Stochastic Shrodinger equation time evolution.
Fields (Attributes)
ntraj::Int
: Number of trajectoriestimes::AbstractVector
: The time list of the evolution.states::Vector{Vector{QuantumObject}}
: The list of result states in each trajectory.expect::Matrix
: The expectation values (averaging all trajectories) corresponding to each time point in times
.expect_all::Array
: The expectation values corresponding to each trajectory and each time point in times
converged::Bool
: Whether the solution is converged or not.alg
: The algorithm which is used during the solving process.abstol::Real
: The absolute tolerance which is used during the solving process.reltol::Real
: The relative tolerance which is used during the solving process.
sourceQuantumToolbox.sesolveProblem
— FunctionsesolveProblem(H::QuantumObject,
ψ0::QuantumObject,
tlist::AbstractVector;
alg::OrdinaryDiffEqAlgorithm=Tsit5()
@@ -300,7 +300,7 @@
H_t::Union{Nothing,Function,TimeDependentOperatorSum}=nothing,
params::NamedTuple=NamedTuple(),
progress_bar::Union{Val,Bool}=Val(true),
- kwargs...)
Generates the ODEProblem for the Schrödinger time evolution of a quantum system:
\[\frac{\partial}{\partial t} |\psi(t)\rangle = -i \hat{H} |\psi(t)\rangle\]
Arguments
H::QuantumObject
: The Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: The initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: The time list of the evolution.alg::OrdinaryDiffEqAlgorithm
: The algorithm used for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: The list of operators to be evaluated during the evolution.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: The time-dependent Hamiltonian of the system. If nothing
, the Hamiltonian is time-independent.params::NamedTuple
: The parameters of the system.progress_bar::Union{Val,Bool}
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.kwargs...
: The keyword arguments passed to the ODEProblem
constructor.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob
: The ODEProblem
for the Schrödinger time evolution of the system.
sourceQuantumToolbox.mesolveProblem
— FunctionmesolveProblem(H::QuantumObject,
+ kwargs...)
Generates the ODEProblem for the Schrödinger time evolution of a quantum system:
\[\frac{\partial}{\partial t} |\psi(t)\rangle = -i \hat{H} |\psi(t)\rangle\]
Arguments
H::QuantumObject
: The Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: The initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: The time list of the evolution.alg::OrdinaryDiffEqAlgorithm
: The algorithm used for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: The list of operators to be evaluated during the evolution.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: The time-dependent Hamiltonian of the system. If nothing
, the Hamiltonian is time-independent.params::NamedTuple
: The parameters of the system.progress_bar::Union{Val,Bool}
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.kwargs...
: The keyword arguments passed to the ODEProblem
constructor.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob
: The ODEProblem
for the Schrödinger time evolution of the system.
sourceQuantumToolbox.mesolveProblem
— FunctionmesolveProblem(H::QuantumObject,
ψ0::QuantumObject,
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
@@ -309,7 +309,7 @@
H_t::Union{Nothing,Function,TimeDependentOperatorSum}=nothing,
params::NamedTuple=NamedTuple(),
progress_bar::Union{Val,Bool}=Val(true),
- kwargs...)
Generates the ODEProblem for the master equation time evolution of an open quantum system:
\[\frac{\partial \hat{\rho}(t)}{\partial t} = -i[\hat{H}, \hat{\rho}(t)] + \sum_n \mathcal{D}(\hat{C}_n) [\hat{\rho}(t)]\]
where
\[\mathcal{D}(\hat{C}_n) [\hat{\rho}(t)] = \hat{C}_n \hat{\rho}(t) \hat{C}_n^\dagger - \frac{1}{2} \hat{C}_n^\dagger \hat{C}_n \hat{\rho}(t) - \frac{1}{2} \hat{\rho}(t) \hat{C}_n^\dagger \hat{C}_n\]
Arguments
H::QuantumObject
: The Hamiltonian $\hat{H}$ or the Liouvillian of the system.ψ0::QuantumObject
: The initial state of the system.tlist::AbstractVector
: The time list of the evolution.c_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm=Tsit5()
: The algorithm used for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the operators for which the expectation values are calculated.H_t::Union{Nothing,Function,TimeDependentOperatorSum}=nothing
: The time-dependent Hamiltonian or Liouvillian.params::NamedTuple=NamedTuple()
: The parameters of the time evolution.progress_bar::Union{Val,Bool}=Val(true)
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.kwargs...
: The keyword arguments for the ODEProblem.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob::ODEProblem
: The ODEProblem for the master equation time evolution.
sourceQuantumToolbox.mcsolveProblem
— FunctionmcsolveProblem(H::QuantumObject{<:AbstractArray{T1},OperatorQuantumObject},
+ kwargs...)
Generates the ODEProblem for the master equation time evolution of an open quantum system:
\[\frac{\partial \hat{\rho}(t)}{\partial t} = -i[\hat{H}, \hat{\rho}(t)] + \sum_n \mathcal{D}(\hat{C}_n) [\hat{\rho}(t)]\]
where
\[\mathcal{D}(\hat{C}_n) [\hat{\rho}(t)] = \hat{C}_n \hat{\rho}(t) \hat{C}_n^\dagger - \frac{1}{2} \hat{C}_n^\dagger \hat{C}_n \hat{\rho}(t) - \frac{1}{2} \hat{\rho}(t) \hat{C}_n^\dagger \hat{C}_n\]
Arguments
H::QuantumObject
: The Hamiltonian $\hat{H}$ or the Liouvillian of the system.ψ0::QuantumObject
: The initial state of the system.tlist::AbstractVector
: The time list of the evolution.c_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm=Tsit5()
: The algorithm used for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the operators for which the expectation values are calculated.H_t::Union{Nothing,Function,TimeDependentOperatorSum}=nothing
: The time-dependent Hamiltonian or Liouvillian.params::NamedTuple=NamedTuple()
: The parameters of the time evolution.progress_bar::Union{Val,Bool}=Val(true)
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.kwargs...
: The keyword arguments for the ODEProblem.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob::ODEProblem
: The ODEProblem for the master equation time evolution.
sourceQuantumToolbox.mcsolveProblem
— FunctionmcsolveProblem(H::QuantumObject{<:AbstractArray{T1},OperatorQuantumObject},
ψ0::QuantumObject{<:AbstractArray{T2},KetQuantumObject},
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
@@ -318,7 +318,7 @@
H_t::Union{Nothing,Function,TimeDependentOperatorSum}=nothing,
params::NamedTuple=NamedTuple(),
jump_callback::TJC=ContinuousLindbladJumpCallback(),
- kwargs...)
Generates the ODEProblem for a single trajectory of the Monte Carlo wave function time evolution of an open quantum system.
Given a system Hamiltonian $\hat{H}$ and a list of collapse (jump) operators $\{\hat{C}_n\}_n$, the evolution of the state $|\psi(t)\rangle$ is governed by the Schrodinger equation:
\[\frac{\partial}{\partial t} |\psi(t)\rangle= -i \hat{H}_{\textrm{eff}} |\psi(t)\rangle\]
with a non-Hermitian effective Hamiltonian:
\[\hat{H}_{\textrm{eff}} = \hat{H} - \frac{i}{2} \sum_n \hat{C}_n^\dagger \hat{C}_n.\]
To the first-order of the wave function in a small time $\delta t$, the strictly negative non-Hermitian portion in $\hat{H}_{\textrm{eff}}$ gives rise to a reduction in the norm of the wave function, namely
\[\langle \psi(t+\delta t) | \psi(t+\delta t) \rangle = 1 - \delta p,\]
where
\[\delta p = \delta t \sum_n \langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle\]
is the corresponding quantum jump probability.
If the environmental measurements register a quantum jump, the wave function undergoes a jump into a state defined by projecting $|\psi(t)\rangle$ using the collapse operator $\hat{C}_n$ corresponding to the measurement, namely
\[| \psi(t+\delta t) \rangle = \frac{\hat{C}_n |\psi(t)\rangle}{ \sqrt{\langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle} }\]
Arguments
H::QuantumObject
: Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: Initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.c_ops::Union{Nothing,AbstractVector,Tuple}
: List of collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.seeds::Union{Nothing, Vector{Int}}
: List of seeds for the random number generator. Length must be equal to the number of trajectories provided.jump_callback::LindbladJumpCallbackType
: The Jump Callback type: Discrete or Continuous.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob::ODEProblem
: The ODEProblem for the Monte Carlo wave function time evolution.
sourceQuantumToolbox.mcsolveEnsembleProblem
— FunctionmcsolveEnsembleProblem(H::QuantumObject{<:AbstractArray{T1},OperatorQuantumObject},
+ kwargs...)
Generates the ODEProblem for a single trajectory of the Monte Carlo wave function time evolution of an open quantum system.
Given a system Hamiltonian $\hat{H}$ and a list of collapse (jump) operators $\{\hat{C}_n\}_n$, the evolution of the state $|\psi(t)\rangle$ is governed by the Schrodinger equation:
\[\frac{\partial}{\partial t} |\psi(t)\rangle= -i \hat{H}_{\textrm{eff}} |\psi(t)\rangle\]
with a non-Hermitian effective Hamiltonian:
\[\hat{H}_{\textrm{eff}} = \hat{H} - \frac{i}{2} \sum_n \hat{C}_n^\dagger \hat{C}_n.\]
To the first-order of the wave function in a small time $\delta t$, the strictly negative non-Hermitian portion in $\hat{H}_{\textrm{eff}}$ gives rise to a reduction in the norm of the wave function, namely
\[\langle \psi(t+\delta t) | \psi(t+\delta t) \rangle = 1 - \delta p,\]
where
\[\delta p = \delta t \sum_n \langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle\]
is the corresponding quantum jump probability.
If the environmental measurements register a quantum jump, the wave function undergoes a jump into a state defined by projecting $|\psi(t)\rangle$ using the collapse operator $\hat{C}_n$ corresponding to the measurement, namely
\[| \psi(t+\delta t) \rangle = \frac{\hat{C}_n |\psi(t)\rangle}{ \sqrt{\langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle} }\]
Arguments
H::QuantumObject
: Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: Initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.c_ops::Union{Nothing,AbstractVector,Tuple}
: List of collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.seeds::Union{Nothing, Vector{Int}}
: List of seeds for the random number generator. Length must be equal to the number of trajectories provided.jump_callback::LindbladJumpCallbackType
: The Jump Callback type: Discrete or Continuous.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob::ODEProblem
: The ODEProblem for the Monte Carlo wave function time evolution.
sourceQuantumToolbox.mcsolveEnsembleProblem
— FunctionmcsolveEnsembleProblem(H::QuantumObject{<:AbstractArray{T1},OperatorQuantumObject},
ψ0::QuantumObject{<:AbstractArray{T2},KetQuantumObject},
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
@@ -329,7 +329,7 @@
jump_callback::TJC=ContinuousLindbladJumpCallback(),
prob_func::Function=_mcsolve_prob_func,
output_func::Function=_mcsolve_output_func,
- kwargs...)
Generates the EnsembleProblem
of ODEProblem
s for the ensemble of trajectories of the Monte Carlo wave function time evolution of an open quantum system.
Given a system Hamiltonian $\hat{H}$ and a list of collapse (jump) operators $\{\hat{C}_n\}_n$, the evolution of the state $|\psi(t)\rangle$ is governed by the Schrodinger equation:
\[\frac{\partial}{\partial t} |\psi(t)\rangle= -i \hat{H}_{\textrm{eff}} |\psi(t)\rangle\]
with a non-Hermitian effective Hamiltonian:
\[\hat{H}_{\textrm{eff}} = \hat{H} - \frac{i}{2} \sum_n \hat{C}_n^\dagger \hat{C}_n.\]
To the first-order of the wave function in a small time $\delta t$, the strictly negative non-Hermitian portion in $\hat{H}_{\textrm{eff}}$ gives rise to a reduction in the norm of the wave function, namely
\[\langle \psi(t+\delta t) | \psi(t+\delta t) \rangle = 1 - \delta p,\]
where
\[\delta p = \delta t \sum_n \langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle\]
is the corresponding quantum jump probability.
If the environmental measurements register a quantum jump, the wave function undergoes a jump into a state defined by projecting $|\psi(t)\rangle$ using the collapse operator $\hat{C}_n$ corresponding to the measurement, namely
\[| \psi(t+\delta t) \rangle = \frac{\hat{C}_n |\psi(t)\rangle}{ \sqrt{\langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle} }\]
Arguments
H::QuantumObject
: Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: Initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.c_ops::Union{Nothing,AbstractVector,Tuple}
: List of collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.seeds::Union{Nothing, Vector{Int}}
: List of seeds for the random number generator. Length must be equal to the number of trajectories provided.jump_callback::LindbladJumpCallbackType
: The Jump Callback type: Discrete or Continuous.prob_func::Function
: Function to use for generating the ODEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob::EnsembleProblem with ODEProblem
: The Ensemble ODEProblem for the Monte Carlo wave function time evolution.
sourceQuantumToolbox.ssesolveProblem
— FunctionssesolveProblem(H::QuantumObject,
+ kwargs...)
Generates the EnsembleProblem
of ODEProblem
s for the ensemble of trajectories of the Monte Carlo wave function time evolution of an open quantum system.
Given a system Hamiltonian $\hat{H}$ and a list of collapse (jump) operators $\{\hat{C}_n\}_n$, the evolution of the state $|\psi(t)\rangle$ is governed by the Schrodinger equation:
\[\frac{\partial}{\partial t} |\psi(t)\rangle= -i \hat{H}_{\textrm{eff}} |\psi(t)\rangle\]
with a non-Hermitian effective Hamiltonian:
\[\hat{H}_{\textrm{eff}} = \hat{H} - \frac{i}{2} \sum_n \hat{C}_n^\dagger \hat{C}_n.\]
To the first-order of the wave function in a small time $\delta t$, the strictly negative non-Hermitian portion in $\hat{H}_{\textrm{eff}}$ gives rise to a reduction in the norm of the wave function, namely
\[\langle \psi(t+\delta t) | \psi(t+\delta t) \rangle = 1 - \delta p,\]
where
\[\delta p = \delta t \sum_n \langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle\]
is the corresponding quantum jump probability.
If the environmental measurements register a quantum jump, the wave function undergoes a jump into a state defined by projecting $|\psi(t)\rangle$ using the collapse operator $\hat{C}_n$ corresponding to the measurement, namely
\[| \psi(t+\delta t) \rangle = \frac{\hat{C}_n |\psi(t)\rangle}{ \sqrt{\langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle} }\]
Arguments
H::QuantumObject
: Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: Initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.c_ops::Union{Nothing,AbstractVector,Tuple}
: List of collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.seeds::Union{Nothing, Vector{Int}}
: List of seeds for the random number generator. Length must be equal to the number of trajectories provided.jump_callback::LindbladJumpCallbackType
: The Jump Callback type: Discrete or Continuous.prob_func::Function
: Function to use for generating the ODEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob::EnsembleProblem with ODEProblem
: The Ensemble ODEProblem for the Monte Carlo wave function time evolution.
sourceQuantumToolbox.ssesolveProblem
— FunctionssesolveProblem(H::QuantumObject,
ψ0::QuantumObject,
tlist::AbstractVector;
sc_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
@@ -341,7 +341,7 @@
```
where
- \]math K = \hat{H} + i \sumn \left(\frac{ej} Cn - \frac{1}{2} \sum{j} Cn^\dagger Cn - \frac{ej^2}{8}\right),
math Mn = Cn - \frac{en}{2}, and
math en = \langle Cn + C_n^\dagger \rangle. ```
Above, C_n
is the n
-th collapse operator and dW_j(t)
is the real Wiener increment associated to C_n
.
Arguments
H::QuantumObject
: The Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: The initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: The time list of the evolution.sc_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: List of stochastic collapse operators $\{\hat{C}_n\}_n$.alg::StochasticDiffEqAlgorithm
: The algorithm used for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of operators to be evaluated during the evolution.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: The time-dependent Hamiltonian of the system. If nothing
, the Hamiltonian is time-independent.params::NamedTuple
: The parameters of the system.kwargs...
: The keyword arguments passed to the SDEProblem
constructor.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-2
and abstol=1e-2
. - For more details about
alg
please refer to DifferentialEquations.jl
(SDE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob
: The SDEProblem
for the Stochastic Schrödinger time evolution of the system.
sourceQuantumToolbox.ssesolveEnsembleProblem
— FunctionssesolveEnsembleProblem(H::QuantumObject,
+ \]math K = \hat{H} + i \sumn \left(\frac{ej} Cn - \frac{1}{2} \sum{j} Cn^\dagger Cn - \frac{ej^2}{8}\right),
math Mn = Cn - \frac{en}{2}, and
math en = \langle Cn + C_n^\dagger \rangle. ```
Above, C_n
is the n
-th collapse operator and dW_j(t)
is the real Wiener increment associated to C_n
.
Arguments
H::QuantumObject
: The Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: The initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: The time list of the evolution.sc_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: List of stochastic collapse operators $\{\hat{C}_n\}_n$.alg::StochasticDiffEqAlgorithm
: The algorithm used for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of operators to be evaluated during the evolution.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: The time-dependent Hamiltonian of the system. If nothing
, the Hamiltonian is time-independent.params::NamedTuple
: The parameters of the system.kwargs...
: The keyword arguments passed to the SDEProblem
constructor.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-2
and abstol=1e-2
. - For more details about
alg
please refer to DifferentialEquations.jl
(SDE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob
: The SDEProblem
for the Stochastic Schrödinger time evolution of the system.
sourceQuantumToolbox.ssesolveEnsembleProblem
— FunctionssesolveEnsembleProblem(H::QuantumObject,
ψ0::QuantumObject,
tlist::AbstractVector;
sc_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
@@ -355,7 +355,7 @@
```
where
- \]math K = \hat{H} + i \sumn \left(\frac{ej} Cn - \frac{1}{2} \sum{j} Cn^\dagger Cn - \frac{ej^2}{8}\right),
math Mn = Cn - \frac{en}{2}, and
math en = \langle Cn + C_n^\dagger \rangle. ```
Above, C_n
is the n
-th collapse operator and dW_j(t)
is the real Wiener increment associated to C_n
.
Arguments
H::QuantumObject
: The Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: The initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: The time list of the evolution.sc_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: List of stochastic collapse operators $\{\hat{C}_n\}_n$.alg::StochasticDiffEqAlgorithm
: The algorithm used for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of operators to be evaluated during the evolution.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: The time-dependent Hamiltonian of the system. If nothing
, the Hamiltonian is time-independent.params::NamedTuple
: The parameters of the system.prob_func::Function
: Function to use for generating the SDEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.kwargs...
: The keyword arguments passed to the SDEProblem
constructor.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-2
and abstol=1e-2
. - For more details about
alg
please refer to DifferentialEquations.jl
(SDE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob::EnsembleProblem with SDEProblem
: The Ensemble SDEProblem for the Stochastic Shrödinger time evolution.
sourceQuantumToolbox.lr_mesolveProblem
— Functionlr_mesolveProblem(H, z, B, tlist, c_ops; e_ops=(), f_ops=(), opt=LRMesolveOptions(), kwargs...) where T
+ \]math K = \hat{H} + i \sumn \left(\frac{ej} Cn - \frac{1}{2} \sum{j} Cn^\dagger Cn - \frac{ej^2}{8}\right),
math Mn = Cn - \frac{en}{2}, and
math en = \langle Cn + C_n^\dagger \rangle. ```
Above, C_n
is the n
-th collapse operator and dW_j(t)
is the real Wiener increment associated to C_n
.
Arguments
H::QuantumObject
: The Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: The initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: The time list of the evolution.sc_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: List of stochastic collapse operators $\{\hat{C}_n\}_n$.alg::StochasticDiffEqAlgorithm
: The algorithm used for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of operators to be evaluated during the evolution.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: The time-dependent Hamiltonian of the system. If nothing
, the Hamiltonian is time-independent.params::NamedTuple
: The parameters of the system.prob_func::Function
: Function to use for generating the SDEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.kwargs...
: The keyword arguments passed to the SDEProblem
constructor.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-2
and abstol=1e-2
. - For more details about
alg
please refer to DifferentialEquations.jl
(SDE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
prob::EnsembleProblem with SDEProblem
: The Ensemble SDEProblem for the Stochastic Shrödinger time evolution.
sourceQuantumToolbox.lr_mesolveProblem
— Functionlr_mesolveProblem(H, z, B, tlist, c_ops; e_ops=(), f_ops=(), opt=LRMesolveOptions(), kwargs...) where T
Formulates the ODEproblem for the low-rank time evolution of the system. The function is called by lr_mesolve.
Parameters
@@ -377,7 +377,7 @@
opt : LRMesolveOptions
The options of the problem.
kwargs : NamedTuple
- Additional keyword arguments for the ODEProblem.
sourceQuantumToolbox.sesolve
— Functionsesolve(H::QuantumObject,
+ Additional keyword arguments for the ODEProblem.
sourceQuantumToolbox.sesolve
— Functionsesolve(H::QuantumObject,
ψ0::QuantumObject,
tlist::AbstractVector;
alg::OrdinaryDiffEqAlgorithm=Tsit5(),
@@ -385,7 +385,7 @@
H_t::Union{Nothing,Function,TimeDependentOperatorSum}=nothing,
params::NamedTuple=NamedTuple(),
progress_bar::Union{Val,Bool}=Val(true),
- kwargs...)
Time evolution of a closed quantum system using the Schrödinger equation:
\[\frac{\partial}{\partial t} |\psi(t)\rangle = -i \hat{H} |\psi(t)\rangle\]
Arguments
H::QuantumObject
: The Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: The initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.progress_bar::Union{Val,Bool}
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionSol
: The solution of the time evolution. See also TimeEvolutionSol
sourceQuantumToolbox.mesolve
— Functionmesolve(H::QuantumObject,
+ kwargs...)
Time evolution of a closed quantum system using the Schrödinger equation:
\[\frac{\partial}{\partial t} |\psi(t)\rangle = -i \hat{H} |\psi(t)\rangle\]
Arguments
H::QuantumObject
: The Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: The initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.progress_bar::Union{Val,Bool}
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionSol
: The solution of the time evolution. See also TimeEvolutionSol
sourceQuantumToolbox.mesolve
— Functionmesolve(H::QuantumObject,
ψ0::QuantumObject,
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
@@ -394,7 +394,7 @@
H_t::Union{Nothing,Function,TimeDependentOperatorSum}=nothing,
params::NamedTuple=NamedTuple(),
progress_bar::Union{Val,Bool}=Val(true),
- kwargs...)
Time evolution of an open quantum system using Lindblad master equation:
\[\frac{\partial \hat{\rho}(t)}{\partial t} = -i[\hat{H}, \hat{\rho}(t)] + \sum_n \mathcal{D}(\hat{C}_n) [\hat{\rho}(t)]\]
where
\[\mathcal{D}(\hat{C}_n) [\hat{\rho}(t)] = \hat{C}_n \hat{\rho}(t) \hat{C}_n^\dagger - \frac{1}{2} \hat{C}_n^\dagger \hat{C}_n \hat{\rho}(t) - \frac{1}{2} \hat{\rho}(t) \hat{C}_n^\dagger \hat{C}_n\]
Arguments
H::QuantumObject
: The Hamiltonian $\hat{H}$ or the Liouvillian of the system.ψ0::QuantumObject
: The initial state of the system.tlist::AbstractVector
: The time list of the evolution.c_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Named Tuple of parameters to pass to the solver.progress_bar::Union{Val,Bool}
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionSol
: The solution of the time evolution. See also TimeEvolutionSol
sourceQuantumToolbox.mcsolve
— Functionmcsolve(H::QuantumObject{<:AbstractArray{T1},OperatorQuantumObject},
+ kwargs...)
Time evolution of an open quantum system using Lindblad master equation:
\[\frac{\partial \hat{\rho}(t)}{\partial t} = -i[\hat{H}, \hat{\rho}(t)] + \sum_n \mathcal{D}(\hat{C}_n) [\hat{\rho}(t)]\]
where
\[\mathcal{D}(\hat{C}_n) [\hat{\rho}(t)] = \hat{C}_n \hat{\rho}(t) \hat{C}_n^\dagger - \frac{1}{2} \hat{C}_n^\dagger \hat{C}_n \hat{\rho}(t) - \frac{1}{2} \hat{\rho}(t) \hat{C}_n^\dagger \hat{C}_n\]
Arguments
H::QuantumObject
: The Hamiltonian $\hat{H}$ or the Liouvillian of the system.ψ0::QuantumObject
: The initial state of the system.tlist::AbstractVector
: The time list of the evolution.c_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Named Tuple of parameters to pass to the solver.progress_bar::Union{Val,Bool}
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionSol
: The solution of the time evolution. See also TimeEvolutionSol
sourceQuantumToolbox.mcsolve
— Functionmcsolve(H::QuantumObject{<:AbstractArray{T1},OperatorQuantumObject},
ψ0::QuantumObject{<:AbstractArray{T2},KetQuantumObject},
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
@@ -410,7 +410,7 @@
output_func::Function = _mcsolve_output_func,
progress_bar::Union{Val,Bool} = Val(true),
kwargs...,
-)
Time evolution of an open quantum system using quantum trajectories.
Given a system Hamiltonian $\hat{H}$ and a list of collapse (jump) operators $\{\hat{C}_n\}_n$, the evolution of the state $|\psi(t)\rangle$ is governed by the Schrodinger equation:
\[\frac{\partial}{\partial t} |\psi(t)\rangle= -i \hat{H}_{\textrm{eff}} |\psi(t)\rangle\]
with a non-Hermitian effective Hamiltonian:
\[\hat{H}_{\textrm{eff}} = \hat{H} - \frac{i}{2} \sum_n \hat{C}_n^\dagger \hat{C}_n.\]
To the first-order of the wave function in a small time $\delta t$, the strictly negative non-Hermitian portion in $\hat{H}_{\textrm{eff}}$ gives rise to a reduction in the norm of the wave function, namely
\[\langle \psi(t+\delta t) | \psi(t+\delta t) \rangle = 1 - \delta p,\]
where
\[\delta p = \delta t \sum_n \langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle\]
is the corresponding quantum jump probability.
If the environmental measurements register a quantum jump, the wave function undergoes a jump into a state defined by projecting $|\psi(t)\rangle$ using the collapse operator $\hat{C}_n$ corresponding to the measurement, namely
\[| \psi(t+\delta t) \rangle = \frac{\hat{C}_n |\psi(t)\rangle}{ \sqrt{\langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle} }\]
Arguments
H::QuantumObject
: Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: Initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.c_ops::Union{Nothing,AbstractVector,Tuple}
: List of collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.seeds::Union{Nothing, Vector{Int}}
: List of seeds for the random number generator. Length must be equal to the number of trajectories provided.ntraj::Int
: Number of trajectories to use.ensemble_method
: Ensemble method to use.jump_callback::LindbladJumpCallbackType
: The Jump Callback type: Discrete or Continuous.prob_func::Function
: Function to use for generating the ODEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.kwargs...
: Additional keyword arguments to pass to the solver.progress_bar::Union{Val,Bool}
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.
Notes
ensemble_method
can be one of EnsembleThreads()
, EnsembleSerial()
, EnsembleDistributed()
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionMCSol
: The solution of the time evolution. See also TimeEvolutionMCSol
sourceQuantumToolbox.ssesolve
— Functionssesolve(H::QuantumObject,
+)
Time evolution of an open quantum system using quantum trajectories.
Given a system Hamiltonian $\hat{H}$ and a list of collapse (jump) operators $\{\hat{C}_n\}_n$, the evolution of the state $|\psi(t)\rangle$ is governed by the Schrodinger equation:
\[\frac{\partial}{\partial t} |\psi(t)\rangle= -i \hat{H}_{\textrm{eff}} |\psi(t)\rangle\]
with a non-Hermitian effective Hamiltonian:
\[\hat{H}_{\textrm{eff}} = \hat{H} - \frac{i}{2} \sum_n \hat{C}_n^\dagger \hat{C}_n.\]
To the first-order of the wave function in a small time $\delta t$, the strictly negative non-Hermitian portion in $\hat{H}_{\textrm{eff}}$ gives rise to a reduction in the norm of the wave function, namely
\[\langle \psi(t+\delta t) | \psi(t+\delta t) \rangle = 1 - \delta p,\]
where
\[\delta p = \delta t \sum_n \langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle\]
is the corresponding quantum jump probability.
If the environmental measurements register a quantum jump, the wave function undergoes a jump into a state defined by projecting $|\psi(t)\rangle$ using the collapse operator $\hat{C}_n$ corresponding to the measurement, namely
\[| \psi(t+\delta t) \rangle = \frac{\hat{C}_n |\psi(t)\rangle}{ \sqrt{\langle \psi(t) | \hat{C}_n^\dagger \hat{C}_n | \psi(t) \rangle} }\]
Arguments
H::QuantumObject
: Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: Initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.c_ops::Union{Nothing,AbstractVector,Tuple}
: List of collapse operators $\{\hat{C}_n\}_n$.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.seeds::Union{Nothing, Vector{Int}}
: List of seeds for the random number generator. Length must be equal to the number of trajectories provided.ntraj::Int
: Number of trajectories to use.ensemble_method
: Ensemble method to use.jump_callback::LindbladJumpCallbackType
: The Jump Callback type: Discrete or Continuous.prob_func::Function
: Function to use for generating the ODEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.kwargs...
: Additional keyword arguments to pass to the solver.progress_bar::Union{Val,Bool}
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.
Notes
ensemble_method
can be one of EnsembleThreads()
, EnsembleSerial()
, EnsembleDistributed()
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-6
and abstol=1e-8
. - For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionMCSol
: The solution of the time evolution. See also TimeEvolutionMCSol
sourceQuantumToolbox.ssesolve
— Functionssesolve(H::QuantumObject,
ψ0::QuantumObject{<:AbstractArray{T2},KetQuantumObject},
tlist::AbstractVector,
sc_ops::Union{Nothing, AbstractVector}=nothing;
@@ -427,7 +427,7 @@
```
where
- \]math K = \hat{H} + i \sumn \left(\frac{ej} Cn - \frac{1}{2} \sum{j} Cn^\dagger Cn - \frac{ej^2}{8}\right),
math Mn = Cn - \frac{en}{2}, and
math en = \langle Cn + C_n^\dagger \rangle. ```
Above, C_n
is the n
-th collapse operator and dW_j(t)
is the real Wiener increment associated to C_n
.
Arguments
H::QuantumObject
: Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: Initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.sc_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: List of stochastic collapse operators $\{\hat{C}_n\}_n$.alg::StochasticDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.seeds::Union{Nothing, Vector{Int}}
: List of seeds for the random number generator. Length must be equal to the number of trajectories provided.ntraj::Int
: Number of trajectories to use.ensemble_method
: Ensemble method to use.prob_func::Function
: Function to use for generating the SDEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.progress_bar::Union{Val,Bool}
: Whether to show a progress bar.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
ensemble_method
can be one of EnsembleThreads()
, EnsembleSerial()
, EnsembleDistributed()
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-2
and abstol=1e-2
. - For more details about
alg
please refer to DifferentialEquations.jl
(SDE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionSSESol
: The solution of the time evolution. See also TimeEvolutionSSESol
sourceQuantumToolbox.dfd_mesolve
— Functiondfd_mesolve(H::Function, ψ0::QuantumObject,
+ \]math K = \hat{H} + i \sumn \left(\frac{ej} Cn - \frac{1}{2} \sum{j} Cn^\dagger Cn - \frac{ej^2}{8}\right),
math Mn = Cn - \frac{en}{2}, and
math en = \langle Cn + C_n^\dagger \rangle. ```
Above, C_n
is the n
-th collapse operator and dW_j(t)
is the real Wiener increment associated to C_n
.
Arguments
H::QuantumObject
: Hamiltonian of the system $\hat{H}$.ψ0::QuantumObject
: Initial state of the system $|\psi(0)\rangle$.tlist::AbstractVector
: List of times at which to save the state of the system.sc_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: List of stochastic collapse operators $\{\hat{C}_n\}_n$.alg::StochasticDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::Union{Nothing,AbstractVector,Tuple}
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function,TimeDependentOperatorSum}
: Time-dependent part of the Hamiltonian.params::NamedTuple
: Dictionary of parameters to pass to the solver.seeds::Union{Nothing, Vector{Int}}
: List of seeds for the random number generator. Length must be equal to the number of trajectories provided.ntraj::Int
: Number of trajectories to use.ensemble_method
: Ensemble method to use.prob_func::Function
: Function to use for generating the SDEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.progress_bar::Union{Val,Bool}
: Whether to show a progress bar.kwargs...
: Additional keyword arguments to pass to the solver.
Notes
ensemble_method
can be one of EnsembleThreads()
, EnsembleSerial()
, EnsembleDistributed()
- The states will be saved depend on the keyword argument
saveat
in kwargs
. - If
e_ops
is specified, the default value of saveat=[tlist[end]]
(only save the final state), otherwise, saveat=tlist
(saving the states corresponding to tlist
). You can also specify e_ops
and saveat
separately. - The default tolerances in
kwargs
are given as reltol=1e-2
and abstol=1e-2
. - For more details about
alg
please refer to DifferentialEquations.jl
(SDE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionSSESol
: The solution of the time evolution. See also TimeEvolutionSSESol
sourceQuantumToolbox.dfd_mesolve
— Functiondfd_mesolve(H::Function, ψ0::QuantumObject,
t_l::AbstractVector, c_ops::Function, maxdims::AbstractVector,
dfd_params::NamedTuple=NamedTuple();
alg::OrdinaryDiffEqAlgorithm=Tsit5(),
@@ -435,7 +435,7 @@
H_t::Union{Nothing,Function,TimeDependentOperatorSum}=nothing,
params::NamedTuple=NamedTuple(),
tol_list::Vector{<:Number}=fill(1e-8, length(maxdims)),
- kwargs...)
Time evolution of an open quantum system using master equation, dynamically changing the dimension of the Hilbert subspaces.
Notes
- For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
sourceQuantumToolbox.dsf_mesolve
— Functiondsf_mesolve(H::Function,
+ kwargs...)
Time evolution of an open quantum system using master equation, dynamically changing the dimension of the Hilbert subspaces.
Notes
- For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
sourceQuantumToolbox.dsf_mesolve
— Functiondsf_mesolve(H::Function,
ψ0::QuantumObject,
t_l::AbstractVector, c_ops::Function,
op_list::Vector{TOl},
@@ -447,7 +447,7 @@
params::NamedTuple=NamedTuple(),
δα_list::Vector{<:Number}=fill(0.2, length(op_list)),
krylov_dim::Int=max(6, min(10, cld(length(ket2dm(ψ0).data), 4))),
- kwargs...)
Time evolution of an open quantum system using master equation and the Dynamical Shifted Fock algorithm.
Notes
- For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
sourceQuantumToolbox.dsf_mcsolve
— Functiondsf_mcsolve(H::Function,
+ kwargs...)
Time evolution of an open quantum system using master equation and the Dynamical Shifted Fock algorithm.
Notes
- For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
sourceQuantumToolbox.dsf_mcsolve
— Functiondsf_mcsolve(H::Function,
ψ0::QuantumObject,
t_l::AbstractVector, c_ops::Function,
op_list::Vector{TOl},
@@ -463,7 +463,7 @@
jump_callback::LindbladJumpCallbackType=ContinuousLindbladJumpCallback(),
krylov_dim::Int=max(6, min(10, cld(length(ket2dm(ψ0).data), 4))),
progress_bar::Union{Bool,Val} = Val(true)
- kwargs...)
Time evolution of a quantum system using the Monte Carlo wave function method and the Dynamical Shifted Fock algorithm.
Notes
- For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
sourceQuantumToolbox.lr_mesolve
— Functionlr_mesolve(prob::ODEProblem; kwargs...)
+ kwargs...)
Time evolution of a quantum system using the Monte Carlo wave function method and the Dynamical Shifted Fock algorithm.
Notes
- For more details about
alg
please refer to DifferentialEquations.jl
(ODE Solvers) - For more details about
kwargs
please refer to DifferentialEquations.jl
(Keyword Arguments)
sourceQuantumToolbox.lr_mesolve
— Functionlr_mesolve(prob::ODEProblem; kwargs...)
Solves the ODEProblem formulated by lr_mesolveProblem. The function is called by lr_mesolve.
Parameters
@@ -471,13 +471,13 @@
prob : ODEProblem
The ODEProblem formulated by lr_mesolveProblem.
kwargs : NamedTuple
- Additional keyword arguments for the ODEProblem.
sourceQuantumToolbox.liouvillian
— Functionliouvillian(H::QuantumObject, c_ops::Union{Nothing,AbstractVector,Tuple}=nothing, Id_cache=I(prod(H.dims)))
Construct the Liouvillian SuperOperator
for a system Hamiltonian $\hat{H}$ and a set of collapse operators $\{\hat{C}_n\}_n$:
\[\mathcal{L} [\cdot] = -i[\hat{H}, \cdot] + \sum_n \mathcal{D}(\hat{C}_n) [\cdot]\]
where
\[\mathcal{D}(\hat{C}_n) [\cdot] = \hat{C}_n [\cdot] \hat{C}_n^\dagger - \frac{1}{2} \hat{C}_n^\dagger \hat{C}_n [\cdot] - \frac{1}{2} [\cdot] \hat{C}_n^\dagger \hat{C}_n\]
The optional argument Id_cache
can be used to pass a precomputed identity matrix. This can be useful when the same function is applied multiple times with a known Hilbert space dimension.
See also spre
, spost
, and lindblad_dissipator
.
sourceQuantumToolbox.liouvillian_generalized
— Functionliouvillian_generalized(H::QuantumObject, fields::Vector,
-T_list::Vector; N_trunc::Int=size(H,1), tol::Float64=0.0, σ_filter::Union{Nothing, Real}=nothing)
Constructs the generalized Liouvillian for a system coupled to a bath of harmonic oscillators.
See, e.g., Settineri, Alessio, et al. "Dissipation and thermal noise in hybrid quantum systems in the ultrastrong-coupling regime." Physical Review A 98.5 (2018): 053834.
sourceQuantumToolbox.steadystate
— Functionsteadystate(
+ Additional keyword arguments for the ODEProblem.
sourceQuantumToolbox.liouvillian
— Functionliouvillian(H::QuantumObject, c_ops::Union{Nothing,AbstractVector,Tuple}=nothing, Id_cache=I(prod(H.dims)))
Construct the Liouvillian SuperOperator
for a system Hamiltonian $\hat{H}$ and a set of collapse operators $\{\hat{C}_n\}_n$:
\[\mathcal{L} [\cdot] = -i[\hat{H}, \cdot] + \sum_n \mathcal{D}(\hat{C}_n) [\cdot]\]
where
\[\mathcal{D}(\hat{C}_n) [\cdot] = \hat{C}_n [\cdot] \hat{C}_n^\dagger - \frac{1}{2} \hat{C}_n^\dagger \hat{C}_n [\cdot] - \frac{1}{2} [\cdot] \hat{C}_n^\dagger \hat{C}_n\]
The optional argument Id_cache
can be used to pass a precomputed identity matrix. This can be useful when the same function is applied multiple times with a known Hilbert space dimension.
See also spre
, spost
, and lindblad_dissipator
.
sourceQuantumToolbox.liouvillian_generalized
— Functionliouvillian_generalized(H::QuantumObject, fields::Vector,
+T_list::Vector; N_trunc::Int=size(H,1), tol::Float64=0.0, σ_filter::Union{Nothing, Real}=nothing)
Constructs the generalized Liouvillian for a system coupled to a bath of harmonic oscillators.
See, e.g., Settineri, Alessio, et al. "Dissipation and thermal noise in hybrid quantum systems in the ultrastrong-coupling regime." Physical Review A 98.5 (2018): 053834.
sourceQuantumToolbox.steadystate
— Functionsteadystate(
H::QuantumObject,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
solver::SteadyStateSolver = SteadyStateDirectSolver(),
kwargs...
-)
Solve the stationary state based on different solvers.
Parameters
H::QuantumObject
: The Hamiltonian or the Liouvillian of the system.c_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the collapse operators.solver::SteadyStateSolver=SteadyStateDirectSolver()
: see documentation Solving for Steady-State Solutions for different solvers.kwargs...
: The keyword arguments for the solver.
sourcesteadystate(
+)
Solve the stationary state based on different solvers.
Parameters
H::QuantumObject
: The Hamiltonian or the Liouvillian of the system.c_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the collapse operators.solver::SteadyStateSolver=SteadyStateDirectSolver()
: see documentation Solving for Steady-State Solutions for different solvers.kwargs...
: The keyword arguments for the solver.
sourcesteadystate(
H::QuantumObject,
ψ0::QuantumObject,
tspan::Real = Inf,
@@ -486,7 +486,7 @@
reltol::Real = 1.0e-8,
abstol::Real = 1.0e-10,
kwargs...
-)
Solve the stationary state based on time evolution (ordinary differential equations; OrdinaryDiffEq.jl
) with a given initial state.
The termination condition of the stationary state $|\rho\rangle\rangle$ is that either the following condition is true
:
\[\lVert\frac{\partial |\hat{\rho}\rangle\rangle}{\partial t}\rVert \leq \textrm{reltol} \times\lVert\frac{\partial |\hat{\rho}\rangle\rangle}{\partial t}+|\hat{\rho}\rangle\rangle\rVert\]
or
\[\lVert\frac{\partial |\hat{\rho}\rangle\rangle}{\partial t}\rVert \leq \textrm{abstol}\]
Parameters
H::QuantumObject
: The Hamiltonian or the Liouvillian of the system.ψ0::QuantumObject
: The initial state of the system.tspan::Real=Inf
: The final time step for the steady state problem.c_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the collapse operators.solver::SteadyStateODESolver=SteadyStateODESolver()
: see SteadyStateODESolver
for more details.reltol::Real=1.0e-8
: Relative tolerance in steady state terminate condition and solver adaptive timestepping.abstol::Real=1.0e-10
: Absolute tolerance in steady state terminate condition and solver adaptive timestepping.kwargs...
: The keyword arguments for the ODEProblem.
sourceQuantumToolbox.steadystate_floquet
— Functionsteadystate_floquet(
+)
Solve the stationary state based on time evolution (ordinary differential equations; OrdinaryDiffEq.jl
) with a given initial state.
The termination condition of the stationary state $|\rho\rangle\rangle$ is that either the following condition is true
:
\[\lVert\frac{\partial |\hat{\rho}\rangle\rangle}{\partial t}\rVert \leq \textrm{reltol} \times\lVert\frac{\partial |\hat{\rho}\rangle\rangle}{\partial t}+|\hat{\rho}\rangle\rangle\rVert\]
or
\[\lVert\frac{\partial |\hat{\rho}\rangle\rangle}{\partial t}\rVert \leq \textrm{abstol}\]
Parameters
H::QuantumObject
: The Hamiltonian or the Liouvillian of the system.ψ0::QuantumObject
: The initial state of the system.tspan::Real=Inf
: The final time step for the steady state problem.c_ops::Union{Nothing,AbstractVector,Tuple}=nothing
: The list of the collapse operators.solver::SteadyStateODESolver=SteadyStateODESolver()
: see SteadyStateODESolver
for more details.reltol::Real=1.0e-8
: Relative tolerance in steady state terminate condition and solver adaptive timestepping.abstol::Real=1.0e-10
: Absolute tolerance in steady state terminate condition and solver adaptive timestepping.kwargs...
: The keyword arguments for the ODEProblem.
sourceQuantumToolbox.steadystate_floquet
— Functionsteadystate_floquet(
H_0::QuantumObject{MT,OpType1},
H_p::QuantumObject{<:AbstractArray,OpType2},
H_m::QuantumObject{<:AbstractArray,OpType3},
@@ -510,7 +510,7 @@
\vdots \\
\hat{\rho}_{n_{\textrm{max}}-1} \\
\hat{\rho}_{n_{\textrm{max}}}
-\end{pmatrix}\]This will allow to simultaneously obtain all the $\hat{\rho}_n$.
In the case of SSFloquetEffectiveLiouvillian
, instead, the effective Liouvillian is calculated using the matrix continued fraction method.
different return The two solvers returns different objects. The SSFloquetLinearSystem
returns a list of QuantumObject
, containing the density matrices for each Fourier component ($\hat{\rho}_{-n}$, with $n$ from $0$ to $n_\textrm{max}$), while the SSFloquetEffectiveLiouvillian
returns only $\hat{\rho}_0$.
Arguments
H_0::QuantumObject
: The Hamiltonian or the Liouvillian of the undriven system.H_p::QuantumObject
: The Hamiltonian or the Liouvillian of the part of the drive that oscillates as $e^{i \omega t}$.H_m::QuantumObject
: The Hamiltonian or the Liouvillian of the part of the drive that oscillates as $e^{-i \omega t}$.ωd::Number
: The frequency of the drive.c_ops::Union{Nothing,AbstractVector} = nothing
: The optional collapse operators.n_max::Integer = 2
: The number of Fourier components to consider.tol::R = 1e-8
: The tolerance for the solver.solver::FSolver = SSFloquetLinearSystem
: The solver to use.kwargs...
: Additional keyword arguments to be passed to the solver.
sourceQuantumToolbox.SteadyStateDirectSolver
— TypeSteadyStateDirectSolver()
A solver which solves steadystate
by finding the inverse of Liouvillian SuperOperator
using the standard method given in LinearAlgebra
.
sourceQuantumToolbox.SteadyStateEigenSolver
— TypeSteadyStateEigenSolver()
A solver which solves steadystate
by finding the zero (or lowest) eigenvalue of Liouvillian SuperOperator
using eigsolve
.
sourceQuantumToolbox.SteadyStateLinearSolver
— TypeSteadyStateLinearSolver(alg = KrylovJL_GMRES(), Pl = nothing, Pr = nothing)
A solver which solves steadystate
by finding the inverse of Liouvillian SuperOperator
using the alg
orithms given in LinearSolve.jl
.
Parameters
alg::SciMLLinearSolveAlgorithm=KrylovJL_GMRES()
: algorithms given in LinearSolve.jl
Pl::Union{Function,Nothing}=nothing
: left preconditioner, see documentation Solving for Steady-State Solutions for more details.Pr::Union{Function,Nothing}=nothing
: right preconditioner, see documentation Solving for Steady-State Solutions for more details.
sourceQuantumToolbox.SteadyStateODESolver
— TypeSteadyStateODESolver(alg = Tsit5())
An ordinary differential equation (ODE) solver for solving steadystate
.
It includes a field (attribute) SteadyStateODESolver.alg
that specifies the solving algorithm. Default to Tsit5()
.
For more details about the solvers, please refer to OrdinaryDiffEq.jl
sourceCorrelations and Spectrum
QuantumToolbox.correlation_3op_2t
— Functioncorrelation_3op_2t(H::QuantumObject,
+\end{pmatrix}\]This will allow to simultaneously obtain all the $\hat{\rho}_n$.
In the case of SSFloquetEffectiveLiouvillian
, instead, the effective Liouvillian is calculated using the matrix continued fraction method.
different return The two solvers returns different objects. The SSFloquetLinearSystem
returns a list of QuantumObject
, containing the density matrices for each Fourier component ($\hat{\rho}_{-n}$, with $n$ from $0$ to $n_\textrm{max}$), while the SSFloquetEffectiveLiouvillian
returns only $\hat{\rho}_0$.
Arguments
H_0::QuantumObject
: The Hamiltonian or the Liouvillian of the undriven system.H_p::QuantumObject
: The Hamiltonian or the Liouvillian of the part of the drive that oscillates as $e^{i \omega t}$.H_m::QuantumObject
: The Hamiltonian or the Liouvillian of the part of the drive that oscillates as $e^{-i \omega t}$.ωd::Number
: The frequency of the drive.c_ops::Union{Nothing,AbstractVector} = nothing
: The optional collapse operators.n_max::Integer = 2
: The number of Fourier components to consider.tol::R = 1e-8
: The tolerance for the solver.solver::FSolver = SSFloquetLinearSystem
: The solver to use.kwargs...
: Additional keyword arguments to be passed to the solver.
sourceQuantumToolbox.SteadyStateDirectSolver
— TypeSteadyStateDirectSolver()
A solver which solves steadystate
by finding the inverse of Liouvillian SuperOperator
using the standard method given in LinearAlgebra
.
sourceQuantumToolbox.SteadyStateEigenSolver
— TypeSteadyStateEigenSolver()
A solver which solves steadystate
by finding the zero (or lowest) eigenvalue of Liouvillian SuperOperator
using eigsolve
.
sourceQuantumToolbox.SteadyStateLinearSolver
— TypeSteadyStateLinearSolver(alg = KrylovJL_GMRES(), Pl = nothing, Pr = nothing)
A solver which solves steadystate
by finding the inverse of Liouvillian SuperOperator
using the alg
orithms given in LinearSolve.jl
.
Parameters
alg::SciMLLinearSolveAlgorithm=KrylovJL_GMRES()
: algorithms given in LinearSolve.jl
Pl::Union{Function,Nothing}=nothing
: left preconditioner, see documentation Solving for Steady-State Solutions for more details.Pr::Union{Function,Nothing}=nothing
: right preconditioner, see documentation Solving for Steady-State Solutions for more details.
sourceQuantumToolbox.SteadyStateODESolver
— TypeSteadyStateODESolver(alg = Tsit5())
An ordinary differential equation (ODE) solver for solving steadystate
.
It includes a field (attribute) SteadyStateODESolver.alg
that specifies the solving algorithm. Default to Tsit5()
.
For more details about the solvers, please refer to OrdinaryDiffEq.jl
sourceCorrelations and Spectrum
QuantumToolbox.correlation_3op_2t
— Functioncorrelation_3op_2t(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector,
τ_l::AbstractVector,
@@ -518,7 +518,7 @@
B::QuantumObject,
C::QuantumObject,
c_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
- kwargs...)
Returns the two-times correlation function of three operators $\hat{A}$, $\hat{B}$ and $\hat{C}$: $\expval{\hat{A}(t) \hat{B}(t + \tau) \hat{C}(t)}$
for a given initial state $\ket{\psi_0}$.
sourceQuantumToolbox.correlation_2op_2t
— Functioncorrelation_2op_2t(H::QuantumObject,
+ kwargs...)
Returns the two-times correlation function of three operators $\hat{A}$, $\hat{B}$ and $\hat{C}$: $\expval{\hat{A}(t) \hat{B}(t + \tau) \hat{C}(t)}$
for a given initial state $\ket{\psi_0}$.
sourceQuantumToolbox.correlation_2op_2t
— Functioncorrelation_2op_2t(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector,
τ_l::AbstractVector,
@@ -526,20 +526,20 @@
B::QuantumObject,
c_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
reverse::Bool=false,
- kwargs...)
Returns the two-times correlation function of two operators $\hat{A}$ and $\hat{B}$ at different times: $\expval{\hat{A}(t + \tau) \hat{B}(t)}$.
When reverse=true
, the correlation function is calculated as $\expval{\hat{A}(t) \hat{B}(t + \tau)}$.
sourceQuantumToolbox.correlation_2op_1t
— Functioncorrelation_2op_1t(H::QuantumObject,
+ kwargs...)
Returns the two-times correlation function of two operators $\hat{A}$ and $\hat{B}$ at different times: $\expval{\hat{A}(t + \tau) \hat{B}(t)}$.
When reverse=true
, the correlation function is calculated as $\expval{\hat{A}(t) \hat{B}(t + \tau)}$.
sourceQuantumToolbox.correlation_2op_1t
— Functioncorrelation_2op_1t(H::QuantumObject,
ψ0::QuantumObject,
τ_l::AbstractVector,
A::QuantumObject,
B::QuantumObject,
c_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
reverse::Bool=false,
- kwargs...)
Returns the one-time correlation function of two operators $\hat{A}$ and $\hat{B}$ at different times $\expval{\hat{A}(\tau) \hat{B}(0)}$.
When reverse=true
, the correlation function is calculated as $\expval{\hat{A}(0) \hat{B}(\tau)}$.
sourceQuantumToolbox.spectrum
— Functionspectrum(H::QuantumObject,
+ kwargs...)
Returns the one-time correlation function of two operators $\hat{A}$ and $\hat{B}$ at different times $\expval{\hat{A}(\tau) \hat{B}(0)}$.
When reverse=true
, the correlation function is calculated as $\expval{\hat{A}(0) \hat{B}(\tau)}$.
sourceQuantumToolbox.spectrum
— Functionspectrum(H::QuantumObject,
ω_list::AbstractVector,
A::QuantumObject{<:AbstractArray{T2},OperatorQuantumObject},
B::QuantumObject{<:AbstractArray{T3},OperatorQuantumObject},
c_ops::Union{Nothing,AbstractVector,Tuple}=nothing;
solver::MySolver=ExponentialSeries(),
- kwargs...)
Returns the emission spectrum
\[S(\omega) = \int_{-\infty}^\infty \expval{\hat{A}(\tau) \hat{B}(0)} e^{-i \omega \tau} d \tau\]
sourceMetrics
QuantumToolbox.entropy_vn
— Functionentropy_vn(ρ::QuantumObject; base::Int=0, tol::Real=1e-15)
Calculates the Von Neumann entropy $S = - \Tr \left[ \hat{\rho} \log \left( \hat{\rho} \right) \right]$ where $\hat{\rho}$ is the density matrix of the system.
The base
parameter specifies the base of the logarithm to use, and when using the default value 0, the natural logarithm is used. The tol
parameter describes the absolute tolerance for detecting the zero-valued eigenvalues of the density matrix $\hat{\rho}$.
Examples
Pure state:
julia> ψ = fock(2,0)
+ kwargs...)
Returns the emission spectrum
\[S(\omega) = \int_{-\infty}^\infty \expval{\hat{A}(\tau) \hat{B}(0)} e^{-i \omega \tau} d \tau\]
sourceMetrics
QuantumToolbox.entropy_vn
— Functionentropy_vn(ρ::QuantumObject; base::Int=0, tol::Real=1e-15)
Calculates the Von Neumann entropy $S = - \Tr \left[ \hat{\rho} \log \left( \hat{\rho} \right) \right]$ where $\hat{\rho}$ is the density matrix of the system.
The base
parameter specifies the base of the logarithm to use, and when using the default value 0, the natural logarithm is used. The tol
parameter describes the absolute tolerance for detecting the zero-valued eigenvalues of the density matrix $\hat{\rho}$.
Examples
Pure state:
julia> ψ = fock(2,0)
Quantum Object: type=Ket dims=[2] size=(2,)
2-element Vector{ComplexF64}:
1.0 + 0.0im
@@ -559,8 +559,8 @@
⋅ 0.5-0.0im
julia> entropy_vn(ρ, base=2)
-1.0
sourceQuantumToolbox.entanglement
— Functionentanglement(QO::QuantumObject, sel::Union{Int,AbstractVector{Int},Tuple})
Calculates the entanglement by doing the partial trace of QO
, selecting only the dimensions with the indices contained in the sel
vector, and then using the Von Neumann entropy entropy_vn
.
sourceQuantumToolbox.tracedist
— Functiontracedist(ρ::QuantumObject, σ::QuantumObject)
Calculates the trace distance between two QuantumObject
: $T(\hat{\rho}, \hat{\sigma}) = \frac{1}{2} \lVert \hat{\rho} - \hat{\sigma} \rVert_1$
sourceQuantumToolbox.fidelity
— Functionfidelity(ρ::QuantumObject, σ::QuantumObject)
Calculate the fidelity of two QuantumObject
: $F(\hat{\rho}, \hat{\sigma}) = \textrm{Tr} \sqrt{\sqrt{\hat{\rho}} \hat{\sigma} \sqrt{\hat{\rho}}}$
Here, the definition is from Nielsen & Chuang, "Quantum Computation and Quantum Information". It is the square root of the fidelity defined in R. Jozsa, Journal of Modern Optics, 41:12, 2315 (1994).
sourceMiscellaneous
QuantumToolbox.wigner
— Functionwigner(state::QuantumObject, xvec::AbstractVector, yvec::AbstractVector; g::Real=√2,
- solver::WignerSolver=WignerLaguerre())
Generates the Wigner quasipropability distribution of state
at points xvec + 1im * yvec
. The g
parameter is a scaling factor related to the value of $\hbar$ in the commutation relation $[x, y] = i \hbar$ via $\hbar=2/g^2$ giving the default value $\hbar=1$.
The solver
parameter can be either WignerLaguerre()
or WignerClenshaw()
. The former uses the Laguerre polynomial expansion of the Wigner function, while the latter uses the Clenshaw algorithm. The Laguerre expansion is faster for sparse matrices, while the Clenshaw algorithm is faster for dense matrices. The WignerLaguerre
solver has an optional parallel
parameter which defaults to true
and uses multithreading to speed up the calculation.
sourceQuantumToolbox.negativity
— Functionnegativity(ρ::QuantumObject, subsys::Int; logarithmic::Bool=false)
Compute the negativity $N(\hat{\rho}) = \frac{\Vert \hat{\rho}^{\Gamma}\Vert_1 - 1}{2}$ where $\hat{\rho}^{\Gamma}$ is the partial transpose of $\hat{\rho}$ with respect to the subsystem, and $\Vert \hat{X} \Vert_1=\textrm{Tr}\sqrt{\hat{X}^\dagger \hat{X}}$ is the trace norm.
Arguments
ρ::QuantumObject
: The density matrix (ρ.type
must be OperatorQuantumObject
).subsys::Int
: an index that indicates which subsystem to compute the negativity for.logarithmic::Bool
: choose whether to calculate logarithmic negativity or not. Default as false
Returns
N::Real
: The value of negativity.
Examples
julia> Ψ = bell_state(0, 0)
+1.0
sourceQuantumToolbox.entanglement
— Functionentanglement(QO::QuantumObject, sel::Union{Int,AbstractVector{Int},Tuple})
Calculates the entanglement by doing the partial trace of QO
, selecting only the dimensions with the indices contained in the sel
vector, and then using the Von Neumann entropy entropy_vn
.
sourceQuantumToolbox.tracedist
— Functiontracedist(ρ::QuantumObject, σ::QuantumObject)
Calculates the trace distance between two QuantumObject
: $T(\hat{\rho}, \hat{\sigma}) = \frac{1}{2} \lVert \hat{\rho} - \hat{\sigma} \rVert_1$
sourceQuantumToolbox.fidelity
— Functionfidelity(ρ::QuantumObject, σ::QuantumObject)
Calculate the fidelity of two QuantumObject
: $F(\hat{\rho}, \hat{\sigma}) = \textrm{Tr} \sqrt{\sqrt{\hat{\rho}} \hat{\sigma} \sqrt{\hat{\rho}}}$
Here, the definition is from Nielsen & Chuang, "Quantum Computation and Quantum Information". It is the square root of the fidelity defined in R. Jozsa, Journal of Modern Optics, 41:12, 2315 (1994).
sourceMiscellaneous
QuantumToolbox.wigner
— Functionwigner(state::QuantumObject, xvec::AbstractVector, yvec::AbstractVector; g::Real=√2,
+ solver::WignerSolver=WignerLaguerre())
Generates the Wigner quasipropability distribution of state
at points xvec + 1im * yvec
. The g
parameter is a scaling factor related to the value of $\hbar$ in the commutation relation $[x, y] = i \hbar$ via $\hbar=2/g^2$ giving the default value $\hbar=1$.
The solver
parameter can be either WignerLaguerre()
or WignerClenshaw()
. The former uses the Laguerre polynomial expansion of the Wigner function, while the latter uses the Clenshaw algorithm. The Laguerre expansion is faster for sparse matrices, while the Clenshaw algorithm is faster for dense matrices. The WignerLaguerre
solver has an optional parallel
parameter which defaults to true
and uses multithreading to speed up the calculation.
sourceQuantumToolbox.negativity
— Functionnegativity(ρ::QuantumObject, subsys::Int; logarithmic::Bool=false)
Compute the negativity $N(\hat{\rho}) = \frac{\Vert \hat{\rho}^{\Gamma}\Vert_1 - 1}{2}$ where $\hat{\rho}^{\Gamma}$ is the partial transpose of $\hat{\rho}$ with respect to the subsystem, and $\Vert \hat{X} \Vert_1=\textrm{Tr}\sqrt{\hat{X}^\dagger \hat{X}}$ is the trace norm.
Arguments
ρ::QuantumObject
: The density matrix (ρ.type
must be OperatorQuantumObject
).subsys::Int
: an index that indicates which subsystem to compute the negativity for.logarithmic::Bool
: choose whether to calculate logarithmic negativity or not. Default as false
Returns
N::Real
: The value of negativity.
Examples
julia> Ψ = bell_state(0, 0)
Quantum Object: type=Ket dims=[2, 2] size=(4,)
4-element Vector{ComplexF64}:
0.7071067811865475 + 0.0im
@@ -577,7 +577,7 @@
0.5+0.0im 0.0+0.0im 0.0+0.0im 0.5+0.0im
julia> negativity(ρ, 2)
-0.4999999999999998
sourceLinear Maps
QuantumToolbox.AbstractLinearMap
— TypeAbstractLinearMap{T, TS}
Represents a general linear map with element type T
and size TS
.
Overview
A linear map is a transformation L
that satisfies:
- Additivity:
math L(u + v) = L(u) + L(v)
- Homogeneity:
math L(cu) = cL(u)
It is typically represented as a matrix with dimensions given by size
, and this abtract type helps to define this map when the matrix is not explicitly available.
Methods
Base.eltype(A)
: Returns the element type T
.Base.size(A)
: Returns the size A.size
.Base.size(A, i)
: Returns the i
-th dimension.
Example
As an example, we now define the linear map used in the eigsolve_al
function for Arnoldi-Lindblad eigenvalue solver:
struct ArnoldiLindbladIntegratorMap{T,TS,TI} <: AbstractLinearMap{T,TS}
+0.4999999999999998
sourceLinear Maps
QuantumToolbox.AbstractLinearMap
— TypeAbstractLinearMap{T, TS}
Represents a general linear map with element type T
and size TS
.
Overview
A linear map is a transformation L
that satisfies:
- Additivity:
math L(u + v) = L(u) + L(v)
- Homogeneity:
math L(cu) = cL(u)
It is typically represented as a matrix with dimensions given by size
, and this abtract type helps to define this map when the matrix is not explicitly available.
Methods
Base.eltype(A)
: Returns the element type T
.Base.size(A)
: Returns the size A.size
.Base.size(A, i)
: Returns the i
-th dimension.
Example
As an example, we now define the linear map used in the eigsolve_al
function for Arnoldi-Lindblad eigenvalue solver:
struct ArnoldiLindbladIntegratorMap{T,TS,TI} <: AbstractLinearMap{T,TS}
elty::Type{T}
size::TS
integrator::TI
@@ -587,15 +587,15 @@
reinit!(A.integrator, x)
solve!(A.integrator)
return copyto!(y, A.integrator.u)
-end
where integrator
is the ODE integrator for the time-evolution. In this way, we can diagonalize this linear map using the eigsolve
function.
sourceUtility functions
QuantumToolbox.versioninfo
— FunctionQuantumToolbox.versioninfo(io::IO=stdout)
Command line output of information on QuantumToolbox, dependencies, and system information, same as QuantumToolbox.about
.
sourceQuantumToolbox.about
— FunctionQuantumToolbox.about(io::IO=stdout)
Command line output of information on QuantumToolbox, dependencies, and system information, same as QuantumToolbox.versioninfo
.
sourceQuantumToolbox.gaussian
— Functiongaussian(x::Number, μ::Number, σ::Number)
Returns the gaussian function $\exp \left[- \frac{(x - \mu)^2}{2 \sigma^2} \right]$, where $\mu$ and $\sigma^2$ are the mean and the variance respectively.
sourceQuantumToolbox.n_thermal
— Functionn_thermal(ω::Real, ω_th::Real)
Return the number of photons in thermal equilibrium for an harmonic oscillator mode with frequency $\omega$, at the temperature described by $\omega_{\textrm{th}} \equiv k_B T / \hbar$:
\[n(\omega, \omega_{\textrm{th}}) = \frac{1}{e^{\omega/\omega_{\textrm{th}}} - 1},\]
where $\hbar$ is the reduced Planck constant, and $k_B$ is the Boltzmann constant.
sourceQuantumToolbox.PhysicalConstants
— Constantconst PhysicalConstants
A NamedTuple
which stores some constant values listed in CODATA recommended values of the fundamental physical constants: 2022
The current stored constants are:
c
: (exact) speed of light in vacuum with unit $[\textrm{m}\cdot\textrm{s}^{-1}]$G
: Newtonian constant of gravitation with unit $[\textrm{m}^3\cdot\textrm{kg}^{−1}\cdot\textrm{s}^{−2}]$h
: (exact) Planck constant with unit $[\textrm{J}\cdot\textrm{s}]$ħ
: reduced Planck constant with unit $[\textrm{J}\cdot\textrm{s}]$e
: (exact) elementary charge with unit $[\textrm{C}]$μ0
: vacuum magnetic permeability with unit $[\textrm{N}\cdot\textrm{A}^{-2}]$ϵ0
: vacuum electric permittivity with unit $[\textrm{F}\cdot\textrm{m}^{-1}]$k
: (exact) Boltzmann constant with unit $[\textrm{J}\cdot\textrm{K}^{-1}]$NA
: (exact) Avogadro constant with unit $[\textrm{mol}^{-1}]$
Examples
julia> PhysicalConstants.ħ
-1.0545718176461565e-34
sourceQuantumToolbox.convert_unit
— Functionconvert_unit(value::Real, unit1::Symbol, unit2::Symbol)
Convert the energy value
from unit1
to unit2
. The unit1
and unit2
can be either the following Symbol
:
:J
: Joule:eV
: electron volt:meV
: milli-electron volt:MHz
: Mega-Hertz multiplied by Planck constant $h$:GHz
: Giga-Hertz multiplied by Planck constant $h$:K
: Kelvin multiplied by Boltzmann constant $k$:mK
: milli-Kelvin multiplied by Boltzmann constant $k$
Note that we use the values stored in PhysicalConstants
to do the conversion.
Examples
julia> convert_unit(1, :eV, :J)
+end
where integrator
is the ODE integrator for the time-evolution. In this way, we can diagonalize this linear map using the eigsolve
function.
sourceUtility functions
QuantumToolbox.versioninfo
— FunctionQuantumToolbox.versioninfo(io::IO=stdout)
Command line output of information on QuantumToolbox, dependencies, and system information, same as QuantumToolbox.about
.
sourceQuantumToolbox.about
— FunctionQuantumToolbox.about(io::IO=stdout)
Command line output of information on QuantumToolbox, dependencies, and system information, same as QuantumToolbox.versioninfo
.
sourceQuantumToolbox.gaussian
— Functiongaussian(x::Number, μ::Number, σ::Number)
Returns the gaussian function $\exp \left[- \frac{(x - \mu)^2}{2 \sigma^2} \right]$, where $\mu$ and $\sigma^2$ are the mean and the variance respectively.
sourceQuantumToolbox.n_thermal
— Functionn_thermal(ω::Real, ω_th::Real)
Return the number of photons in thermal equilibrium for an harmonic oscillator mode with frequency $\omega$, at the temperature described by $\omega_{\textrm{th}} \equiv k_B T / \hbar$:
\[n(\omega, \omega_{\textrm{th}}) = \frac{1}{e^{\omega/\omega_{\textrm{th}}} - 1},\]
where $\hbar$ is the reduced Planck constant, and $k_B$ is the Boltzmann constant.
sourceQuantumToolbox.PhysicalConstants
— Constantconst PhysicalConstants
A NamedTuple
which stores some constant values listed in CODATA recommended values of the fundamental physical constants: 2022
The current stored constants are:
c
: (exact) speed of light in vacuum with unit $[\textrm{m}\cdot\textrm{s}^{-1}]$G
: Newtonian constant of gravitation with unit $[\textrm{m}^3\cdot\textrm{kg}^{−1}\cdot\textrm{s}^{−2}]$h
: (exact) Planck constant with unit $[\textrm{J}\cdot\textrm{s}]$ħ
: reduced Planck constant with unit $[\textrm{J}\cdot\textrm{s}]$e
: (exact) elementary charge with unit $[\textrm{C}]$μ0
: vacuum magnetic permeability with unit $[\textrm{N}\cdot\textrm{A}^{-2}]$ϵ0
: vacuum electric permittivity with unit $[\textrm{F}\cdot\textrm{m}^{-1}]$k
: (exact) Boltzmann constant with unit $[\textrm{J}\cdot\textrm{K}^{-1}]$NA
: (exact) Avogadro constant with unit $[\textrm{mol}^{-1}]$
Examples
julia> PhysicalConstants.ħ
+1.0545718176461565e-34
sourceQuantumToolbox.convert_unit
— Functionconvert_unit(value::Real, unit1::Symbol, unit2::Symbol)
Convert the energy value
from unit1
to unit2
. The unit1
and unit2
can be either the following Symbol
:
:J
: Joule:eV
: electron volt:meV
: milli-electron volt:MHz
: Mega-Hertz multiplied by Planck constant $h$:GHz
: Giga-Hertz multiplied by Planck constant $h$:K
: Kelvin multiplied by Boltzmann constant $k$:mK
: milli-Kelvin multiplied by Boltzmann constant $k$
Note that we use the values stored in PhysicalConstants
to do the conversion.
Examples
julia> convert_unit(1, :eV, :J)
1.602176634e-19
julia> convert_unit(1, :GHz, :J)
6.62607015e-25
julia> convert_unit(1, :meV, :mK)
-11604.518121550082
sourceQuantumToolbox.row_major_reshape
— Functionrow_major_reshape(Q::AbstractArray, shapes...)
Reshapes Q
in the row-major order, as numpy.
sourceQuantumToolbox.meshgrid
— Functionmeshgrid(x::AbstractVector, y::AbstractVector)
Equivalent to numpy meshgrid.
sourceQuantumToolbox._calculate_expectation!
— Function_calculate_expectation!(p,z,B,idx) where T
+11604.518121550082
sourceQuantumToolbox.row_major_reshape
— Functionrow_major_reshape(Q::AbstractArray, shapes...)
Reshapes Q
in the row-major order, as numpy.
sourceQuantumToolbox.meshgrid
— Functionmeshgrid(x::AbstractVector, y::AbstractVector)
Equivalent to numpy meshgrid.
sourceQuantumToolbox._calculate_expectation!
— Function_calculate_expectation!(p,z,B,idx) where T
Calculates the expectation values and function values of the operators and functions in p.e_ops and p.f_ops, respectively, and stores them in p.expvals and p.funvals.
The function is called by the callback _save_affect_lr_mesolve!.
@@ -608,7 +608,7 @@
B : AbstractMatrix{T}
The B matrix.
idx : Integer
- The index of the current time step.
sourceQuantumToolbox._adjM_condition_variational
— Function_adjM_condition_variational(u, t, integrator) where T
+ The index of the current time step.
sourceQuantumToolbox._adjM_condition_variational
— Function_adjM_condition_variational(u, t, integrator) where T
Condition for the dynamical rank adjustment based on the leakage out of the low-rank manifold.
Parameters
@@ -618,14 +618,14 @@
t : Real
The current time.
integrator : ODEIntegrator
- The integrator of the problem.
sourceQuantumToolbox._adjM_affect!
— Function_adjM_affect!(integrator)
+ The integrator of the problem.
sourceQuantumToolbox._adjM_affect!
— Function_adjM_affect!(integrator)
Affect function for the dynamical rank adjustment. It increases the rank of the low-rank manifold by one, and updates the matrices accordingly.
If Δt>0, it rewinds the integrator to the previous time step.
Parameters
----------
integrator : ODEIntegrator
- The integrator of the problem.
sourceQuantumToolbox._adjM_condition_ratio
— Function_adjM_condition_ratio(u, t, integrator) where T
+ The integrator of the problem.
sourceQuantumToolbox._adjM_condition_ratio
— Function_adjM_condition_ratio(u, t, integrator) where T
Condition for the dynamical rank adjustment based on the ratio between the smallest and largest eigenvalues of the density matrix.
The spectrum of the density matrix is calculated efficiently using the properties of the SVD decomposition of the matrix.
@@ -636,7 +636,7 @@
t : Real
The current time.
integrator : ODEIntegrator
- The integrator of the problem.
sourceQuantumToolbox._pinv!
— Function_pinv!(A, T1, T2; atol::Real=0.0, rtol::Real=(eps(real(float(oneunit(T))))*min(size(A)...))*iszero(atol)) where T
+ The integrator of the problem.
sourceQuantumToolbox._pinv!
— Function_pinv!(A, T1, T2; atol::Real=0.0, rtol::Real=(eps(real(float(oneunit(T))))*min(size(A)...))*iszero(atol)) where T
Computes the pseudo-inverse of a matrix A, and stores it in T1. If T2 is provided, it is used as a temporary matrix.
The algorithm is based on the SVD decomposition of A, and is taken from the Julia package LinearAlgebra.
The difference with respect to the original function is that the cutoff is done with a smooth function instead of a step function.
@@ -651,7 +651,7 @@
atol : Real
Absolute tolerance for the calculation of the pseudo-inverse.
rtol : Real
- Relative tolerance for the calculation of the pseudo-inverse.
sourceQuantumToolbox.dBdz!
— FunctiondBdz!(du, u, p, t) where T
+ Relative tolerance for the calculation of the pseudo-inverse.
sourceQuantumToolbox.dBdz!
— FunctiondBdz!(du, u, p, t) where T
Dynamical evolution equations for the low-rank manifold. The function is called by the ODEProblem.
Parameters
@@ -663,4 +663,4 @@
p : NamedTuple
The parameters of the problem.
t : Real
- The current time.
sourceSettings
This document was generated with Documenter.jl version 1.7.0 on Saturday 5 October 2024. Using Julia version 1.10.5.
+ The current time.