API
Table of contents
Quantum object (Qobj) and type
QuantumToolbox.Space Type
struct Space <: AbstractSpace
size::Int
end
A structure that describes a single Hilbert space with size = size
.
QuantumToolbox.Dimensions Type
struct Dimensions{N,T<:Tuple} <: AbstractDimensions{N}
to::T
end
A structure that describes the Hilbert Space
of each subsystems.
QuantumToolbox.GeneralDimensions Type
struct GeneralDimensions{N,T1<:Tuple,T2<:Tuple} <: AbstractDimensions{N}
to::T1
from::T2
end
A structure that describes the left-hand side (to
) and right-hand side (from
) Hilbert Space
of an Operator
.
QuantumToolbox.AbstractQuantumObject Type
abstract type AbstractQuantumObject{ObjType,DimType,DataType}
Abstract type for all quantum objects like QuantumObject
and QuantumObjectEvolution
.
Example
julia> sigmax() isa AbstractQuantumObject
true
QuantumToolbox.BraQuantumObject Type
BraQuantumObject <: QuantumObjectType
Constructor representing a bra state
QuantumToolbox.Bra Constant
const Bra = BraQuantumObject()
A constant representing the type of BraQuantumObject
: a bra state
QuantumToolbox.KetQuantumObject Type
KetQuantumObject <: QuantumObjectType
Constructor representing a ket state
QuantumToolbox.Ket Constant
const Ket = KetQuantumObject()
A constant representing the type of KetQuantumObject
: a ket state
QuantumToolbox.OperatorQuantumObject Type
OperatorQuantumObject <: QuantumObjectType
Constructor representing an operator
QuantumToolbox.Operator Constant
const Operator = OperatorQuantumObject()
A constant representing the type of OperatorQuantumObject
: an operator
QuantumToolbox.OperatorBraQuantumObject Type
OperatorBraQuantumObject <: QuantumObjectType
Constructor representing a bra state in the SuperOperator
formalism
QuantumToolbox.OperatorBra Constant
const OperatorBra = OperatorBraQuantumObject()
A constant representing the type of OperatorBraQuantumObject
: a bra state in the SuperOperator
formalism
QuantumToolbox.OperatorKetQuantumObject Type
OperatorKetQuantumObject <: QuantumObjectType
Constructor representing a ket state in the SuperOperator
formalism
QuantumToolbox.OperatorKet Constant
const OperatorKet = OperatorKetQuantumObject()
A constant representing the type of OperatorKetQuantumObject
: a ket state in the SuperOperator
formalism
QuantumToolbox.SuperOperatorQuantumObject Type
SuperOperatorQuantumObject <: QuantumObjectType
Constructor representing a super-operator
QuantumToolbox.SuperOperator Constant
const SuperOperator = SuperOperatorQuantumObject()
A constant representing the type of SuperOperatorQuantumObject
: a super-operator
QuantumToolbox.QuantumObject Type
struct QuantumObject{ObjType<:QuantumObjectType,DimType<:AbstractDimensions,DataType<:AbstractArray} <: AbstractQuantumObject{ObjType,DimType,DataType}
data::DataType
type::ObjType
dimensions::DimType
end
Julia structure representing any time-independent quantum objects. For time-dependent cases, see QuantumObjectEvolution
.
dims
property
For a given H::QuantumObject
, H.dims
or getproperty(H, :dims)
returns its dimensions
in the type of integer-vector.
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:
⎡⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠈⠢⡀⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠈⠢⡀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠈⠢⡀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢⎦
julia> a isa QuantumObject
true
julia> a.dims
1-element SVector{1, Int64} with indices SOneTo(1):
20
julia> a.dimensions
Dimensions{1, Tuple{Space}}((Space(20),))
QuantumToolbox.QuantumObjectEvolution Type
struct QuantumObjectEvolution{ObjType<:QuantumObjectType,DimType<:AbstractDimensions,DataType<:AbstractSciMLOperator} <: AbstractQuantumObject{ObjType,DimType,DataType}
data::DataType
type::ObjType
dimensions::DimType
end
Julia struct representing any time-dependent quantum object. The data
field is a AbstractSciMLOperator
object that represents the time-dependent quantum object. It can be seen as
where p
and time t
, and
For time-independent cases, see QuantumObject
, and for more information about AbstractSciMLOperator
, see the SciML documentation.
dims
property
For a given H::QuantumObjectEvolution
, H.dims
or getproperty(H, :dims)
returns its dimensions
in the type of integer-vector.
Examples
This operator can be initialized in the same way as the QuTiP QobjEvo
object. For example
julia> a = tensor(destroy(10), qeye(2))
Quantum Object: type=Operator dims=[10, 2] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 18 stored entries:
⎡⠀⠑⢄⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠀⠑⢄⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠑⢄⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠀⠑⢄⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠀⠑⎦
julia> coef1(p, t) = exp(-1im * t)
coef1 (generic function with 1 method)
julia> op = QobjEvo(a, coef1)
Quantum Object Evo.: type=Operator dims=[10, 2] size=(20, 20) ishermitian=true isconstant=false
ScalarOperator(0.0 + 0.0im) * MatrixOperator(20 × 20)
If there are more than 2 operators, we need to put each set of operator and coefficient function into a two-element Tuple
, and put all these Tuple
s together in a larger Tuple
:
julia> σm = tensor(qeye(10), sigmam())
Quantum Object: type=Operator dims=[10, 2] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 10 stored entries:
⎡⠂⡀⠀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠂⡀⠀⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠂⡀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠂⡀⠀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠂⡀⎦
julia> coef2(p, t) = sin(t)
coef2 (generic function with 1 method)
julia> op1 = QobjEvo(((a, coef1), (σm, coef2)))
Quantum Object Evo.: type=Operator dims=[10, 2] size=(20, 20) ishermitian=true isconstant=false
(ScalarOperator(0.0 + 0.0im) * MatrixOperator(20 × 20) + ScalarOperator(0.0 + 0.0im) * MatrixOperator(20 × 20))
We can also concretize the operator at a specific time t
julia> op1(0.1)
Quantum Object: type=Operator dims=[10, 2] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 28 stored entries:
⎡⠂⡑⢄⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠂⡑⢄⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠂⡑⢄⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠂⡑⢄⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠂⡑⎦
It also supports parameter-dependent time evolution
julia> coef1(p, t) = exp(-1im * p.ω1 * t)
coef1 (generic function with 1 method)
julia> coef2(p, t) = sin(p.ω2 * t)
coef2 (generic function with 1 method)
julia> op1 = QobjEvo(((a, coef1), (σm, coef2)))
Quantum Object Evo.: type=Operator dims=[10, 2] size=(20, 20) ishermitian=true isconstant=false
(ScalarOperator(0.0 + 0.0im) * MatrixOperator(20 × 20) + ScalarOperator(0.0 + 0.0im) * MatrixOperator(20 × 20))
julia> p = (ω1 = 1.0, ω2 = 0.5)
(ω1 = 1.0, ω2 = 0.5)
julia> op1(p, 0.1)
Quantum Object: type=Operator dims=[10, 2] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 28 stored entries:
⎡⠂⡑⢄⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠂⡑⢄⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠂⡑⢄⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠂⡑⢄⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠂⡑⎦
Base.size Function
size(A::AbstractQuantumObject)
size(A::AbstractQuantumObject, idx::Int)
shape(A::AbstractQuantumObject)
shape(A::AbstractQuantumObject, idx::Int)
Returns a tuple containing each dimensions of the array in the AbstractQuantumObject
.
Optionally, you can specify an index (idx
) to just get the corresponding dimension of the array.
Note
shape
is a synonym of size
.
Base.eltype Function
eltype(A::AbstractQuantumObject)
Returns the elements type of the matrix or vector corresponding to the AbstractQuantumObject
A
.
Base.length Function
length(A::AbstractQuantumObject)
Returns the length of the matrix or vector corresponding to the AbstractQuantumObject
A
.
SciMLOperators.cache_operator Function
SciMLOperators.cached_operator(L::AbstractQuantumObject, u)
Allocate caches for AbstractQuantumObject
L
for in-place evaluation with u
-like input vectors.
Here, u
can be in either the following types:
AbstractVector
Ket
-typeQuantumObject
(ifL
is anOperator
)OperatorKet
-typeQuantumObject
(ifL
is aSuperOperator
)
Qobj boolean functions
QuantumToolbox.isbra Function
isbra(A)
Checks if the QuantumObject
A
is a BraQuantumObject
. Default case returns false
for any other inputs.
QuantumToolbox.isket Function
isket(A)
Checks if the QuantumObject
A
is a KetQuantumObject
. Default case returns false
for any other inputs.
QuantumToolbox.isoper Function
isoper(A)
Checks if the AbstractQuantumObject
A
is a OperatorQuantumObject
. Default case returns false
for any other inputs.
QuantumToolbox.isoperbra Function
isoperbra(A)
Checks if the QuantumObject
A
is a OperatorBraQuantumObject
. Default case returns false
for any other inputs.
QuantumToolbox.isoperket Function
isoperket(A)
Checks if the QuantumObject
A
is a OperatorKetQuantumObject
. Default case returns false
for any other inputs.
QuantumToolbox.issuper Function
issuper(A)
Checks if the AbstractQuantumObject
A
is a SuperOperatorQuantumObject
. Default case returns false
for any other inputs.
LinearAlgebra.ishermitian Function
ishermitian(A::AbstractQuantumObject)
isherm(A::AbstractQuantumObject)
Test whether the AbstractQuantumObject
is Hermitian.
Note
isherm
is a synonym of ishermitian
.
LinearAlgebra.issymmetric Function
issymmetric(A::AbstractQuantumObject)
Test whether the AbstractQuantumObject
is symmetric.
LinearAlgebra.isposdef Function
isposdef(A::AbstractQuantumObject)
Test whether the AbstractQuantumObject
is positive definite (and Hermitian) by trying to perform a Cholesky factorization of A
.
QuantumToolbox.isunitary Function
isunitary(U::QuantumObject; kwargs...)
Test whether the QuantumObject
Base.isapprox
to test whether
Note that all the keyword arguments will be passed to Base.isapprox
.
SciMLOperators.iscached Function
SciMLOperators.iscached(A::AbstractQuantumObject)
Test whether the AbstractQuantumObject
A
has preallocated caches for inplace evaluations.
SciMLOperators.isconstant Function
SciMLOperators.isconstant(A::AbstractQuantumObject)
Test whether the AbstractQuantumObject
A
is constant in time. For a QuantumObject
, this function returns true
, while for a QuantumObjectEvolution
, this function returns true
if the operator is constant in time.
Qobj arithmetic and attributes
Base.zero Function
zero(A::AbstractQuantumObject)
Return a similar AbstractQuantumObject
with dims
and type
are same as A
, but data
is a zero-array.
Base.one Function
one(A::AbstractQuantumObject)
Return a similar AbstractQuantumObject
with dims
and type
are same as A
, but data
is an identity matrix.
Note that A
must be Operator
or SuperOperator
.
Base.conj Function
conj(A::AbstractQuantumObject)
Return the element-wise complex conjugation of the AbstractQuantumObject
.
Base.transpose Function
transpose(A::AbstractQuantumObject)
Lazy matrix transpose of the AbstractQuantumObject
.
Base.adjoint Function
A'
adjoint(A::AbstractQuantumObject)
dag(A::AbstractQuantumObject)
Lazy adjoint (conjugate transposition) of the AbstractQuantumObject
Note
A'
and dag(A)
are synonyms of adjoint(A)
.
LinearAlgebra.dot Function
A ⋅ B
dot(A::QuantumObject, B::QuantumObject)
Compute the dot product between two QuantumObject
:
Note that A
and B
should be Ket
or OperatorKet
Note
A ⋅ B
(where ⋅
can be typed by tab-completing \cdot
in the REPL) is a synonym of dot(A, B)
.
dot(i::QuantumObject, A::AbstractQuantumObject j::QuantumObject)
matrix_element(i::QuantumObject, A::AbstractQuantumObject j::QuantumObject)
Compute the generalized dot product dot(i, A*j)
between a AbstractQuantumObject
and two QuantumObject
(i
and j
), namely
Supports the following inputs:
A
is in the type ofSuperOperator
, withi
andj
are bothOperatorKet
Note
matrix_element(i, A, j)
is a synonym of dot(i, A, j)
.
Base.sqrt Function
√(A)
sqrt(A::QuantumObject)
Matrix square root of QuantumObject
Note
√(A)
(where √
can be typed by tab-completing \sqrt
in the REPL) is a synonym of sqrt(A)
.
Base.log Function
log(A::QuantumObject)
Matrix logarithm of QuantumObject
Note that this function only supports for Operator
and SuperOperator
Base.exp Function
exp(A::QuantumObject)
Matrix exponential of QuantumObject
Note that this function only supports for Operator
and SuperOperator
Base.sin Function
sin(A::QuantumObject)
Matrix sine of QuantumObject
, defined as
Note that this function only supports for Operator
and SuperOperator
Base.cos Function
cos(A::QuantumObject)
Matrix cosine of QuantumObject
, defined as
Note that this function only supports for Operator
and SuperOperator
LinearAlgebra.tr Function
tr(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:
⎡⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠈⠢⡀⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠈⠢⡀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠈⠢⡀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢⎦
julia> tr(a' * a)
190.0 + 0.0im
LinearAlgebra.svdvals Function
svdvals(A::QuantumObject)
Return the singular values of a QuantumObject
in descending order
LinearAlgebra.norm Function
norm(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 eitherKet
,Bra
,OperatorKet
, orOperatorBra
, returns the standard vectorp
-norm (defaultp=2
) ofA
.If
A
is eitherOperator
orSuperOperator
, returns Schattenp
-norm (defaultp=1
) ofA
.
Examples
julia> ψ = fock(10, 2)
Quantum Object: type=Ket dims=[10] size=(10,)
10-element Vector{ComplexF64}:
0.0 + 0.0im
0.0 + 0.0im
1.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
julia> norm(ψ)
1.0
LinearAlgebra.normalize Function
normalize(A::QuantumObject, p::Real)
unit(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
unit
is a synonym of normalize
.
Also, see norm
about its definition for different types of QuantumObject
.
LinearAlgebra.normalize! Function
normalize!(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
.
Base.inv Function
inv(A::AbstractQuantumObject)
Matrix inverse of the AbstractQuantumObject
. If A
is a QuantumObjectEvolution
, the inverse is computed at the last computed time.
LinearAlgebra.diag Function
diag(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
QuantumToolbox.proj Function
proj(ψ::QuantumObject)
Return the projector for a Ket
or Bra
type of QuantumObject
QuantumToolbox.ptrace Function
ptrace(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
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
1.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
julia> ptrace(ψ, 2)
Quantum Object: type=Operator dims=[2] size=(2, 2) ishermitian=true
2×2 Matrix{ComplexF64}:
0.0+0.0im 0.0+0.0im
0.0+0.0im 1.0+0.0im
or in an entangled state
julia> ψ = 1 / √2 * (kron(fock(2,0), fock(2,0)) + kron(fock(2,1), fock(2,1)))
Quantum Object: type=Ket dims=[2, 2] size=(4,)
4-element Vector{ComplexF64}:
0.7071067811865475 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.7071067811865475 + 0.0im
julia> ptrace(ψ, 1)
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
QuantumToolbox.purity Function
purity(ρ::QuantumObject)
Calculate the purity of a QuantumObject
:
Note that this function only supports for Ket
, Bra
, and Operator
SparseArrays.permute Function
permute(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:
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.
QuantumToolbox.tidyup Function
tidyup(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
.
QuantumToolbox.tidyup! Function
tidyup!(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
.
QuantumToolbox.get_data Function
get_data(A::AbstractQuantumObject)
Returns the data of a AbstractQuantumObject
.
QuantumToolbox.get_coherence Function
get_coherence(ψ::QuantumObject)
Get the coherence value
It returns both
QuantumToolbox.partial_transpose Function
partial_transpose(ρ::QuantumObject, mask::Vector{Bool})
Return the partial transpose of a density matrix 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 beOperatorQuantumObject
).mask::Vector{Bool}
: A boolean vector selects which subsystems should be transposed.
Returns
ρ_pt::QuantumObject
: The density matrix with the selected subsystems transposed.
Qobj eigenvalues and eigenvectors
QuantumToolbox.EigsolveResult Type
struct EigsolveResult
A struct containing the eigenvalues, the eigenvectors, and some information from the solver
Fields (Attributes)
values::AbstractVector
: the eigenvaluesvectors::AbstractMatrix
: the transformation matrix (eigenvectors)type::Union{Nothing,QuantumObjectType}
: the type ofQuantumObject
, ornothing
means solving eigen equation for general matrixdimensions::Union{Nothing,AbstractDimensions}
: thedimensions
ofQuantumObject
, ornothing
means solving eigen equation for general matrixiter::Int
: the number of iteration during the solving processnumops::Int
: number of times the linear map was applied in krylov methodsconverged::Bool
: Whether the result is converged
dims
property
For a given res::EigsolveResult
, res.dims
or getproperty(res, :dims)
returns its dimensions
in the type of integer-vector.
Examples
One can obtain the eigenvalues and the corresponding QuantumObject
-type eigenvectors by:
julia> result = eigenstates(sigmax())
EigsolveResult: type=Operator dims=[2]
values:
2-element Vector{ComplexF64}:
-1.0 + 0.0im
1.0 + 0.0im
vectors:
2×2 Matrix{ComplexF64}:
-0.707107+0.0im 0.707107+0.0im
0.707107+0.0im 0.707107+0.0im
julia> λ, ψ, U = result;
julia> λ
2-element Vector{ComplexF64}:
-1.0 + 0.0im
1.0 + 0.0im
julia> ψ
2-element Vector{QuantumObject{KetQuantumObject, Dimensions{1, Tuple{Space}}, Vector{ComplexF64}}}:
Quantum Object: type=Ket dims=[2] size=(2,)
2-element Vector{ComplexF64}:
-0.7071067811865475 + 0.0im
0.7071067811865475 + 0.0im
Quantum Object: type=Ket dims=[2] size=(2,)
2-element Vector{ComplexF64}:
0.7071067811865475 + 0.0im
0.7071067811865475 + 0.0im
julia> U
2×2 Matrix{ComplexF64}:
-0.707107+0.0im 0.707107+0.0im
0.707107+0.0im 0.707107+0.0im
QuantumToolbox.eigenenergies Function
eigenenergies(A::QuantumObject; sparse::Bool=false, kwargs...)
Calculate the eigenenergies
Arguments
A::QuantumObject
: theQuantumObject
to solve eigenvaluessparse::Bool
: iffalse
calleigvals(A::QuantumObject; kwargs...)
, otherwise calleigsolve
. Default tofalse
.kwargs
: Additional keyword arguments passed to the solver. Ifsparse=true
, the keyword arguments are passed toeigsolve
, otherwise toeigen
.
Returns
::Vector{<:Number}
: a list of eigenvalues
QuantumToolbox.eigenstates Function
eigenstates(A::QuantumObject; sparse::Bool=false, kwargs...)
Calculate the eigenvalues and corresponding eigenvectors
Arguments
A::QuantumObject
: theQuantumObject
to solve eigenvalues and eigenvectorssparse::Bool
: iffalse
calleigen(A::QuantumObject; kwargs...)
, otherwise calleigsolve
. Default tofalse
.kwargs
: Additional keyword arguments passed to the solver. Ifsparse=true
, the keyword arguments are passed toeigsolve
, otherwise toeigen
.
Returns
::EigsolveResult
: containing the eigenvalues, the eigenvectors, and some information from the solver. see alsoEigsolveResult
LinearAlgebra.eigen Function
LinearAlgebra.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';
julia> E, ψ, U = eigen(H)
EigsolveResult: type=Operator dims=[5]
values:
5-element Vector{ComplexF64}:
-2.8569700138728 + 0.0im
-1.3556261799742608 + 0.0im
1.3322676295501878e-15 + 0.0im
1.3556261799742677 + 0.0im
2.8569700138728056 + 0.0im
vectors:
5×5 Matrix{ComplexF64}:
0.106101+0.0im -0.471249-0.0im … 0.471249-0.0im 0.106101-0.0im
-0.303127-0.0im 0.638838+0.0im 0.638838+0.0im 0.303127-0.0im
0.537348+0.0im -0.279149-0.0im 0.279149-0.0im 0.537348-0.0im
-0.638838-0.0im -0.303127-0.0im -0.303127-0.0im 0.638838+0.0im
0.447214+0.0im 0.447214+0.0im -0.447214-0.0im 0.447214-0.0im
julia> expect(H, ψ[1]) ≈ E[1]
true
LinearAlgebra.eigvals Function
LinearAlgebra.eigvals(A::QuantumObject; kwargs...)
Same as eigen(A::QuantumObject; kwargs...)
but for only the eigenvalues.
QuantumToolbox.eigsolve Function
eigsolve(A::QuantumObject;
v0::Union{Nothing,AbstractVector}=nothing,
sigma::Union{Nothing, Real}=nothing,
eigvals::Int = 1,
krylovdim::Int = max(20, 2*k+1),
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.
Arguments
A::QuantumObject
: theQuantumObject
to solve eigenvalues and eigenvectors.v0::Union{Nothing,AbstractVector}
: the initial vector for the Arnoldi method. Default is a random vector.sigma::Union{Nothing, Real}
: the shift for the eigenvalue problem. Default isnothing
.eigvals::Int
: the number of eigenvalues to compute. Default is1
.krylovdim::Int
: the dimension of the Krylov subspace. Default ismax(20, 2*k+1)
.tol::Real
: the tolerance for the Arnoldi method. Default is1e-8
.maxiter::Int
: the maximum number of iterations for the Arnoldi method. Default is200
.solver::Union{Nothing, SciMLLinearSolveAlgorithm}
: the linear solver algorithm. Default isnothing
.kwargs
: Additional keyword arguments passed to the solver.
Notes
- For more details about
solver
and extrakwargs
, please refer toLinearSolve.jl
Returns
EigsolveResult
: A struct containing the eigenvalues, the eigenvectors, and some information about the eigsolver
QuantumToolbox.eigsolve_al Function
eigsolve_al(
H::Union{AbstractQuantumObject{HOpType},Tuple},
T::Real,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
alg::OrdinaryDiffEqAlgorithm = Tsit5(),
params::NamedTuple = NamedTuple(),
ρ0::AbstractMatrix = rand_dm(prod(H.dimensions)).data,
eigvals::Int = 1,
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. It can be aQuantumObject
, aQuantumObjectEvolution
, or a tuple of the form supported bymesolve
.T
: The time at which to evaluate the time evolution.c_ops
: A vector of collapse operators. Default isnothing
meaning the system is closed.alg
: The differential equation solver algorithm. Default isTsit5()
.params
: ANamedTuple
containing the parameters of the system.ρ0
: The initial density matrix. If not specified, a random density matrix is used.eigvals
: The number of eigenvalues to compute.krylovdim
: The dimension of the Krylov subspace.maxiter
: The maximum number of iterations for the eigsolver.eigstol
: The tolerance for the eigsolver.kwargs
: Additional keyword arguments passed to the differential equation solver.
Notes
For more details about
alg
please refer toDifferentialEquations.jl
(ODE Solvers)For more details about
kwargs
please refer toDifferentialEquations.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.
Qobj manipulation
QuantumToolbox.ket2dm Function
ket2dm(ψ::QuantumObject)
Transform the ket state
QuantumToolbox.expect Function
expect(O::Union{AbstractQuantumObject,Vector{AbstractQuantumObject}}, ψ::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
If ψ
is a density matrix (Operator
), the function calculates
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)
.
List of observables and states
The observable O
and state ψ
can be given as a list of QuantumObject
, it returns a list of expectation values. If both of them are given as a list, it returns a Matrix
of expectation values.
Examples
julia> ψ1 = 1 / √2 * (fock(10,2) + fock(10,4));
julia> ψ2 = coherent(10, 0.6 + 0.8im);
julia> a = destroy(10);
julia> expect(a' * a, ψ1) |> round
3.0 + 0.0im
julia> expect(Hermitian(a' * a), ψ1) |> round
3.0
julia> round.(expect([a' * a, a' + a, a], [ψ1, ψ2]), digits = 1)
3×2 Matrix{ComplexF64}:
3.0+0.0im 1.0+0.0im
0.0+0.0im 1.2-0.0im
0.0+0.0im 0.6+0.8im
QuantumToolbox.variance Function
variance(O::QuantumObject, ψ::Union{QuantumObject,Vector{QuantumObject}})
Variance of the Operator
O
:
where 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.
Base.kron Function
kron(A::AbstractQuantumObject, B::AbstractQuantumObject, ...)
tensor(A::AbstractQuantumObject, B::AbstractQuantumObject, ...)
⊗(A::AbstractQuantumObject, B::AbstractQuantumObject, ...)
A ⊗ B
Returns the Kronecker product
Note
tensor
and ⊗
(where ⊗
can be typed by tab-completing \otimes
in the REPL) are synonyms of kron
.
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:
⎡⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠈⠢⡀⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠈⠢⡀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠈⠢⡀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢⎦
julia> O = kron(a, a);
julia> size(a), size(O)
((20, 20), (400, 400))
julia> a.dims, O.dims
([20], [20, 20])
QuantumToolbox.to_dense Function
to_dense(A::QuantumObject)
Converts a sparse QuantumObject to a dense QuantumObject.
QuantumToolbox.to_sparse Function
to_sparse(A::QuantumObject)
Converts a dense QuantumObject to a sparse QuantumObject.
QuantumToolbox.vec2mat Function
vec2mat(A::AbstractVector)
Converts a vector to a matrix.
vec2mat(A::QuantumObject)
vector_to_operator(A::QuantumObject)
Convert a quantum object from vector (OperatorKetQuantumObject
-type) to matrix (OperatorQuantumObject
-type)
Note
vector_to_operator
is a synonym of vec2mat
.
QuantumToolbox.mat2vec Function
mat2vec(A::QuantumObject)
operator_to_vector(A::QuantumObject)
Convert a quantum object from matrix (OperatorQuantumObject
-type) to vector (OperatorKetQuantumObject
-type)
Note
operator_to_vector
is a synonym of mat2vec
.
mat2vec(A::AbstractMatrix)
Converts a matrix to a vector.
Generate states and operators
QuantumToolbox.zero_ket Function
zero_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{Dimensions,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.
QuantumToolbox.fock Function
fock(N::Int, j::Int=0; dims::Union{Int,AbstractVector{Int},Tuple}=N, sparse::Union{Bool,Val}=Val(false))
Generates a fock state 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.
QuantumToolbox.basis Function
basis(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.
QuantumToolbox.coherent Function
coherent(N::Int, α::Number)
Generates a coherent state
This state is constructed via the displacement operator displace
and zero-fock state fock
:
QuantumToolbox.rand_ket Function
rand_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{Dimensions,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.
QuantumToolbox.fock_dm Function
fock_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.
QuantumToolbox.coherent_dm Function
coherent_dm(N::Int, α::Number)
Density matrix representation of a coherent state.
Constructed via outer product of coherent
.
QuantumToolbox.thermal_dm Function
thermal_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}
: Iftrue
, 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.
QuantumToolbox.maximally_mixed_dm Function
maximally_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{Dimensions,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.
QuantumToolbox.rand_dm Function
rand_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{Dimensions,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
QuantumToolbox.spin_state Function
spin_state(j::Real, m::Real)
Generate the spin state:
The eigenstate of the Spin-j
m
, where where j
is the spin quantum number and can be a non-negative integer or half-integer
See also jmat
.
QuantumToolbox.spin_coherent Function
spin_coherent(j::Real, θ::Real, ϕ::Real)
Generate the coherent spin state (rotation of the
where the rotation operator is defined as
and 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
QuantumToolbox.bell_state Function
bell_state(x::Union{Int}, z::Union{Int})
Return the Bell state depending on the arguments (x, z)
:
(0, 0)
:(0, 1)
:(1, 0)
:(1, 1)
:
Here, x = 1
(z = 1
) means applying Pauli-
Example
julia> bell_state(0, 0)
Quantum Object: type=Ket dims=[2, 2] size=(4,)
4-element Vector{ComplexF64}:
0.7071067811865475 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.7071067811865475 + 0.0im
julia> bell_state(Val(1), Val(0))
Quantum Object: type=Ket dims=[2, 2] size=(4,)
4-element Vector{ComplexF64}:
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.
QuantumToolbox.triplet_states Function
triplet_states()
Return a list of the two particle triplet states:
QuantumToolbox.w_state Function
w_state(n::Union{Int,Val})
Returns the n
-qubit W-state:
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.
QuantumToolbox.ghz_state Function
ghz_state(n::Union{Int,Val}; d::Int=2)
Returns the generalized n
-qudit Greenberger–Horne–Zeilinger (GHZ) state:
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.
QuantumToolbox.rand_unitary Function
rand_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{Dimensions,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, where 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.
QuantumToolbox.jmat Function
jmat(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
::y
::z
::+
::-
:
Note that if the parameter which
is not specified, returns a set of Spin-j
operators:
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
0.5+0.0im ⋅
julia> jmat(0.5, Val(:-))
Quantum Object: type=Operator dims=[2] size=(2, 2) ishermitian=false
2×2 SparseMatrixCSC{ComplexF64, Int64} with 1 stored entry:
⋅ ⋅
1.0+0.0im ⋅
julia> jmat(1.5, Val(:z))
Quantum Object: type=Operator dims=[4] size=(4, 4) ishermitian=true
4×4 SparseMatrixCSC{ComplexF64, Int64} with 4 stored entries:
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.
QuantumToolbox.spin_Jx Function
spin_Jx(j::Real)
j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
QuantumToolbox.spin_Jy Function
spin_Jy(j::Real)
j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
QuantumToolbox.spin_Jz Function
spin_Jz(j::Real)
j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
QuantumToolbox.spin_Jm Function
spin_Jm(j::Real)
j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
QuantumToolbox.spin_Jp Function
spin_Jp(j::Real)
j
, where j
is the spin quantum number and can be a non-negative integer or half-integer.
See also jmat
.
QuantumToolbox.spin_J_set Function
spin_J_set(j::Real)
A set of Spin-j
operators 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
.
QuantumToolbox.destroy Function
destroy(N::Int)
Bosonic annihilation operator with Hilbert space cutoff N
.
This operator acts on a fock state as
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:
⎡⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠈⠢⡀⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠈⠢⡀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠈⠢⡀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢⎦
julia> fock(20, 3)' * a * fock(20, 4)
2.0 + 0.0im
QuantumToolbox.create Function
create(N::Int)
Bosonic creation operator with Hilbert space cutoff N
.
This operator acts on a fock state as
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:
⎡⠢⡀⠀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠈⠢⡀⠀⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠈⠢⡀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠈⠢⡀⠀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠈⠢⡀⎦
julia> fock(20, 4)' * a_d * fock(20, 3)
2.0 + 0.0im
QuantumToolbox.displace Function
displace(N::Int, α::Number)
Generate a displacement operator:
where
QuantumToolbox.squeeze Function
squeeze(N::Int, z::Number)
Generate a single-mode squeeze operator:
where
QuantumToolbox.num Function
num(N::Int)
Bosonic number operator with Hilbert space cutoff N
.
This operator is defined as
QuantumToolbox.position Function
position(N::Int)
Position operator with Hilbert space cutoff N
.
This operator is defined as
QuantumToolbox.momentum Function
momentum(N::Int)
Momentum operator with Hilbert space cutoff N
.
This operator is defined as
QuantumToolbox.phase Function
phase(N::Int, ϕ0::Real=0)
Single-mode Pegg-Barnett phase operator with Hilbert space cutoff
This operator is defined as
where
and
Reference
QuantumToolbox.fdestroy Function
fdestroy(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
The site index j
should satisfy: 1 ≤ j ≤ N
.
Note that we put
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.
QuantumToolbox.fcreate Function
fcreate(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
The site index j
should satisfy: 1 ≤ j ≤ N
.
Note that we put
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.
QuantumToolbox.tunneling Function
tunneling(N::Int, m::Int=1; sparse::Union{Bool,Val{<:Bool}}=Val(false))
Generate a tunneling operator defined as:
where
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.
QuantumToolbox.qft Function
qft(dimensions)
Generates a discrete Fourier transform matrix dimensions
.
The dimensions
can be either the following types:
dimensions::Int
: Number of basis states in the Hilbert space.dimensions::Union{Dimensions,AbstractVector{Int},Tuple}
: list of dimensions representing the each number of basis in the subsystems.
where
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.
QuantumToolbox.eye Function
eye(N::Int; type=Operator, dims=nothing)
qeye(N::Int; type=Operator, dims=nothing)
Identity operator 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
Note
qeye
is a synonym of eye
.
QuantumToolbox.projection Function
projection(N::Int, i::Int, j::Int)
Generates the projection operator N
.
QuantumToolbox.commutator Function
commutator(A::QuantumObject, B::QuantumObject; anti::Bool=false)
Return the commutator (or anti
-commutator) of the two QuantumObject
:
commutator (
anti=false
):anticommutator (
anti=true
):
Note that A
and B
must be Operator
QuantumToolbox.spre Function
spre(A::AbstractQuantumObject, Id_cache=I(size(A,1)))
Returns the SuperOperator
form of A
acting on the left of the density matrix operator:
Since the density matrix is vectorized in OperatorKet
form: SuperOperator
is always a matrix
(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.
QuantumToolbox.spost Function
spost(B::AbstractQuantumObject, Id_cache=I(size(B,1)))
Returns the SuperOperator
form of B
acting on the right of the density matrix operator:
Since the density matrix is vectorized in OperatorKet
form: SuperOperator
is always a matrix
(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.
QuantumToolbox.sprepost Function
sprepost(A::AbstractQuantumObject, B::AbstractQuantumObject)
Returns the SuperOperator
form of A
and B
acting on the left and right of the density matrix operator, respectively:
Since the density matrix is vectorized in OperatorKet
form: SuperOperator
is always a matrix
(see the section in documentation: Superoperators and Vectorized Operators for more details)
QuantumToolbox.lindblad_dissipator Function
lindblad_dissipator(O::AbstractQuantumObject, Id_cache=I(size(O,1))
Returns the Lindblad SuperOperator
defined as
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.
Synonyms of functions for Qobj
QuantumToolbox.Qobj Type
Qobj(A; kwargs...)
Note
Qobj
is a synonym for generating QuantumObject
. See the docstring of QuantumObject
for more details.
QuantumToolbox.QobjEvo Type
QobjEvo(args...; kwargs...)
Note
QobjEvo
is a synonym for generating QuantumObjectEvolution
. See the docstrings of QuantumObjectEvolution
for more details.
QuantumToolbox.shape Function
size(A::AbstractQuantumObject)
size(A::AbstractQuantumObject, idx::Int)
shape(A::AbstractQuantumObject)
shape(A::AbstractQuantumObject, idx::Int)
Returns a tuple containing each dimensions of the array in the AbstractQuantumObject
.
Optionally, you can specify an index (idx
) to just get the corresponding dimension of the array.
Note
shape
is a synonym of size
.
QuantumToolbox.isherm Function
ishermitian(A::AbstractQuantumObject)
isherm(A::AbstractQuantumObject)
Test whether the AbstractQuantumObject
is Hermitian.
Note
isherm
is a synonym of ishermitian
.
QuantumToolbox.trans Function
transpose(A::AbstractQuantumObject)
Lazy matrix transpose of the AbstractQuantumObject
.
QuantumToolbox.dag Function
A'
adjoint(A::AbstractQuantumObject)
dag(A::AbstractQuantumObject)
Lazy adjoint (conjugate transposition) of the AbstractQuantumObject
Note
A'
and dag(A)
are synonyms of adjoint(A)
.
QuantumToolbox.matrix_element Function
matrix_element(i::QuantumObject, A::QuantumObject j::QuantumObject)
Compute the generalized dot product dot(i, A*j)
between three QuantumObject
:
Note that this function is same as dot(i, A, j)
Supports the following inputs:
A
is in the type ofSuperOperator
, withi
andj
are bothOperatorKet
QuantumToolbox.unit Function
normalize(A::QuantumObject, p::Real)
unit(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
unit
is a synonym of normalize
.
Also, see norm
about its definition for different types of QuantumObject
.
QuantumToolbox.tensor Function
kron(A::AbstractQuantumObject, B::AbstractQuantumObject, ...)
tensor(A::AbstractQuantumObject, B::AbstractQuantumObject, ...)
⊗(A::AbstractQuantumObject, B::AbstractQuantumObject, ...)
A ⊗ B
Returns the Kronecker product
Note
tensor
and ⊗
(where ⊗
can be typed by tab-completing \otimes
in the REPL) are synonyms of kron
.
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:
⎡⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠈⠢⡀⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠈⠢⡀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠈⠢⡀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢⎦
julia> O = kron(a, a);
julia> size(a), size(O)
((20, 20), (400, 400))
julia> a.dims, O.dims
([20], [20, 20])
QuantumToolbox.:⊗ Function
kron(A::AbstractQuantumObject, B::AbstractQuantumObject, ...)
tensor(A::AbstractQuantumObject, B::AbstractQuantumObject, ...)
⊗(A::AbstractQuantumObject, B::AbstractQuantumObject, ...)
A ⊗ B
Returns the Kronecker product
Note
tensor
and ⊗
(where ⊗
can be typed by tab-completing \otimes
in the REPL) are synonyms of kron
.
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:
⎡⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀⎤
⎢⠀⠀⠈⠢⡀⠀⠀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠈⠢⡀⠀⠀⠀⎥
⎢⠀⠀⠀⠀⠀⠀⠈⠢⡀⠀⎥
⎣⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢⎦
julia> O = kron(a, a);
julia> size(a), size(O)
((20, 20), (400, 400))
julia> a.dims, O.dims
([20], [20, 20])
QuantumToolbox.qeye Function
eye(N::Int; type=Operator, dims=nothing)
qeye(N::Int; type=Operator, dims=nothing)
Identity operator 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
Note
qeye
is a synonym of eye
.
QuantumToolbox.vector_to_operator Function
vec2mat(A::AbstractVector)
Converts a vector to a matrix.
vec2mat(A::QuantumObject)
vector_to_operator(A::QuantumObject)
Convert a quantum object from vector (OperatorKetQuantumObject
-type) to matrix (OperatorQuantumObject
-type)
Note
vector_to_operator
is a synonym of vec2mat
.
QuantumToolbox.operator_to_vector Function
mat2vec(A::QuantumObject)
operator_to_vector(A::QuantumObject)
Convert a quantum object from matrix (OperatorQuantumObject
-type) to vector (OperatorKetQuantumObject
-type)
Note
operator_to_vector
is a synonym of mat2vec
.
mat2vec(A::AbstractMatrix)
Converts a matrix to a vector.
QuantumToolbox.sqrtm Function
sqrtm(A::QuantumObject)
Matrix square root of Operator
type of QuantumObject
Note that for other types of QuantumObject
use sprt(A)
instead.
QuantumToolbox.logm Function
logm(A::QuantumObject)
Matrix logarithm of QuantumObject
Note that this function is same as log(A)
and only supports for Operator
and SuperOperator
.
QuantumToolbox.expm Function
expm(A::QuantumObject)
Matrix exponential of QuantumObject
Note that this function is same as exp(A)
and only supports for Operator
and SuperOperator
.
QuantumToolbox.sinm Function
sinm(A::QuantumObject)
Matrix sine of QuantumObject
, defined as
Note that this function is same as sin(A)
and only supports for Operator
and SuperOperator
.
QuantumToolbox.cosm Function
cosm(A::QuantumObject)
Matrix cosine of QuantumObject
, defined as
Note that this function is same as cos(A)
and only supports for Operator
and SuperOperator
.
Time evolution
QuantumToolbox.TimeEvolutionProblem Type
struct TimeEvolutionProblem
A Julia constructor for handling the ODEProblem
of the time evolution of quantum systems.
Fields (Attributes)
prob::AbstractSciMLProblem
: TheODEProblem
of the time evolution.times::Abstractvector
: The time list of the evolution.dimensions::AbstractDimensions
: The dimensions of the Hilbert space.kwargs::KWT
: Generic keyword arguments.
dims
property
For a given prob::TimeEvolutionProblem
, prob.dims
or getproperty(prob, :dims)
returns its dimensions
in the type of integer-vector.
QuantumToolbox.TimeEvolutionSol Type
struct 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::Union{AbstractMatrix,Nothing}
: The expectation values corresponding to each time point intimes
.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.
QuantumToolbox.TimeEvolutionMCSol Type
struct 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::Union{AbstractMatrix,Nothing}
: The expectation values (averaging all trajectories) corresponding to each time point intimes
.average_expect::Union{AbstractMatrix,Nothing}
: The expectation values (averaging all trajectories) corresponding to each time point intimes
.runs_expect::Union{AbstractArray,Nothing}
: The expectation values corresponding to each trajectory and each time point intimes
col_times::Vector{Vector{Real}}
: The time records of every quantum jump occurred in each trajectory.col_which::Vector{Vector{Int}}
: The indices of which collapse operator was responsible for each quantum jump incol_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.
QuantumToolbox.TimeEvolutionStochasticSol Type
struct TimeEvolutionStochasticSol
A structure storing the results and some information from solving trajectories of the Stochastic 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::Union{AbstractMatrix,Nothing}
: The expectation values (averaging all trajectories) corresponding to each time point intimes
.average_expect::Union{AbstractMatrix,Nothing}
: The expectation values (averaging all trajectories) corresponding to each time point intimes
.runs_expect::Union{AbstractArray,Nothing}
: The expectation values corresponding to each trajectory and each time point intimes
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.
QuantumToolbox.sesolveProblem Function
sesolveProblem(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject{KetQuantumObject},
tlist::AbstractVector;
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
progress_bar::Union{Val,Bool} = Val(true),
inplace::Union{Val,Bool} = Val(true),
kwargs...,
)
Generate the ODEProblem for the Schrödinger time evolution of a quantum system:
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. tlist
: List of times at which to save either the state or the expectation values of the system.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
: Parameters to pass to the solver. This argument is usually expressed as aNamedTuple
orAbstractVector
of parameters. For more advanced usage, any custom struct can be used.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.inplace
: Whether to use the inplace version of the ODEProblem. The default isVal(true)
. It is recommended to useVal(true)
for better performance, but it is sometimes necessary to useVal(false)
, for example when performing automatic differentiation using Zygote.jl.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-6
andabstol=1e-8
.For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Returns
prob
: TheTimeEvolutionProblem
containing theODEProblem
for the Schrödinger time evolution of the system.
QuantumToolbox.mesolveProblem Function
mesolveProblem(
H::Union{AbstractQuantumObject,Tuple},
ψ0::QuantumObject,
tlist,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
progress_bar::Union{Val,Bool} = Val(true),
inplace::Union{Val,Bool} = Val(true),
kwargs...,
)
Generate the ODEProblem for the master equation time evolution of an open quantum system:
where
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. It can be either a Ket
or aOperator
.tlist
: List of times at which to save either the state or the expectation values of the system.c_ops
: List of collapse operators. It can be either a Vector
or aTuple
.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
: Parameters to pass to the solver. This argument is usually expressed as aNamedTuple
orAbstractVector
of parameters. For more advanced usage, any custom struct can be used.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.inplace
: Whether to use the inplace version of the ODEProblem. The default isVal(true)
. It is recommended to useVal(true)
for better performance, but it is sometimes necessary to useVal(false)
, for example when performing automatic differentiation using Zygote.jl.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-6
andabstol=1e-8
.For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Returns
prob::ODEProblem
: The ODEProblem for the master equation time evolution.
QuantumToolbox.mcsolveProblem Function
mcsolveProblem(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject{KetQuantumObject},
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
rng::AbstractRNG = default_rng(),
jump_callback::TJC = ContinuousLindbladJumpCallback(),
kwargs...,
)
Generate the ODEProblem for a single trajectory of the Monte Carlo wave function time evolution of an open quantum system.
Given a system Hamiltonian
with a non-Hermitian effective Hamiltonian:
To the first-order of the wave function in a small time
where
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
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. tlist
: List of times at which to save either the state or the expectation values of the system.c_ops
: List of collapse operators. It can be either a Vector
or aTuple
.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
: Parameters to pass to the solver. This argument is usually expressed as aNamedTuple
orAbstractVector
of parameters. For more advanced usage, any custom struct can be used.rng
: Random number generator for reproducibility.jump_callback
: The Jump Callback type: Discrete or Continuous. The default isContinuousLindbladJumpCallback()
, which is more precise.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-6
andabstol=1e-8
.For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Returns
prob
: TheTimeEvolutionProblem
containing theODEProblem
for the Monte Carlo wave function time evolution.
QuantumToolbox.mcsolveEnsembleProblem Function
mcsolveEnsembleProblem(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject{KetQuantumObject},
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
rng::AbstractRNG = default_rng(),
ntraj::Int = 500,
ensemblealg::EnsembleAlgorithm = EnsembleThreads(),
jump_callback::TJC = ContinuousLindbladJumpCallback(),
progress_bar::Union{Val,Bool} = Val(true),
prob_func::Union{Function, Nothing} = nothing,
output_func::Union{Tuple,Nothing} = nothing,
kwargs...,
)
Generate 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
with a non-Hermitian effective Hamiltonian:
To the first-order of the wave function in a small time
where
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
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. tlist
: List of times at which to save either the state or the expectation values of the system.c_ops
: List of collapse operators. It can be either a Vector
or aTuple
.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
: Parameters to pass to the solver. This argument is usually expressed as aNamedTuple
orAbstractVector
of parameters. For more advanced usage, any custom struct can be used.rng
: Random number generator for reproducibility.ntraj
: Number of trajectories to use.ensemblealg
: Ensemble algorithm to use. Default toEnsembleThreads()
.jump_callback
: The Jump Callback type: Discrete or Continuous. The default isContinuousLindbladJumpCallback()
, which is more precise.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.prob_func
: Function to use for generating the ODEProblem.output_func
: aTuple
containing theFunction
to use for generating the output of a single trajectory, the (optional)ProgressBar
object, and the (optional)RemoteChannel
object.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-6
andabstol=1e-8
.For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Returns
prob
: TheTimeEvolutionProblem
containing the EnsembleODEProblem
for the Monte Carlo wave function time evolution.
QuantumToolbox.ssesolveProblem Function
ssesolveProblem(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject{KetQuantumObject},
tlist::AbstractVector,
sc_ops::Union{Nothing,AbstractVector,Tuple,AbstractQuantumObject} = nothing;
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
rng::AbstractRNG = default_rng(),
progress_bar::Union{Val,Bool} = Val(true),
store_measurement::Union{Val, Bool} = Val(false),
kwargs...,
)
Generate the SDEProblem for the Stochastic Schrödinger time evolution of a quantum system. This is defined by the following stochastic differential equation:
where
and
Above,
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. tlist
: List of times at which to save either the state or the expectation values of the system.sc_ops
: List of stochastic collapse operators. It can be either a Vector
, aTuple
or aAbstractQuantumObject
. It is recommended to use the last case when only one operator is provided.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
:NullParameters
of parameters to pass to the solver.rng
: Random number generator for reproducibility.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.store_measurement
: Whether to store the measurement results. Default isVal(false)
.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-2
andabstol=1e-2
.For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Performance Tip
When sc_ops
contains only a single operator, it is recommended to pass only that operator as the argument. This ensures that the stochastic noise is diagonal, making the simulation faster.
Returns
prob
: TheSDEProblem
for the Stochastic Schrödinger time evolution of the system.
QuantumToolbox.ssesolveEnsembleProblem Function
ssesolveEnsembleProblem(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject{KetQuantumObject},
tlist::AbstractVector,
sc_ops::Union{Nothing,AbstractVector,Tuple,AbstractQuantumObject} = nothing;
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
rng::AbstractRNG = default_rng(),
ntraj::Int = 500,
ensemblealg::EnsembleAlgorithm = EnsembleThreads(),
prob_func::Union{Function, Nothing} = nothing,
output_func::Union{Tuple,Nothing} = nothing,
progress_bar::Union{Val,Bool} = Val(true),
store_measurement::Union{Val,Bool} = Val(false),
kwargs...,
)
Generate the SDE EnsembleProblem for the Stochastic Schrödinger time evolution of a quantum system. This is defined by the following stochastic differential equation:
where
and
Above,
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. tlist
: List of times at which to save either the state or the expectation values of the system.sc_ops
: List of stochastic collapse operators. It can be either a Vector
, aTuple
or aAbstractQuantumObject
. It is recommended to use the last case when only one operator is provided.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
:NullParameters
of parameters to pass to the solver.rng
: Random number generator for reproducibility.ntraj
: Number of trajectories to use. Default is500
.ensemblealg
: Ensemble method to use. Default toEnsembleThreads()
.jump_callback
: The Jump Callback type: Discrete or Continuous. The default isContinuousLindbladJumpCallback()
, which is more precise.prob_func
: Function to use for generating the SDEProblem.output_func
: aTuple
containing theFunction
to use for generating the output of a single trajectory, the (optional)ProgressBar
object, and the (optional)RemoteChannel
object.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.store_measurement
: Whether to store the measurement results. Default isVal(false)
.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-2
andabstol=1e-2
.For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Performance Tip
When sc_ops
contains only a single operator, it is recommended to pass only that operator as the argument. This ensures that the stochastic noise is diagonal, making the simulation faster.
Returns
prob::EnsembleProblem with SDEProblem
: The Ensemble SDEProblem for the Stochastic Shrödinger time evolution.
QuantumToolbox.smesolveProblem Function
smesolveProblem(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject,
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
sc_ops::Union{Nothing,AbstractVector,Tuple,AbstractQuantumObject} = nothing;
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
rng::AbstractRNG = default_rng(),
progress_bar::Union{Val,Bool} = Val(true),
store_measurement::Union{Val, Bool} = Val(false),
kwargs...,
)
Generate the SDEProblem for the Stochastic Master Equation time evolution of an open quantum system. This is defined by the following stochastic differential equation:
where
is the Lindblad superoperator, and
Above,
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. It can be either a Ket
or aOperator
.tlist
: List of times at which to save either the state or the expectation values of the system.c_ops
: List of collapse operators. It can be either a Vector
or aTuple
.sc_ops
: List of stochastic collapse operators. It can be either a Vector
, aTuple
or aAbstractQuantumObject
. It is recommended to use the last case when only one operator is provided.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
:NullParameters
of parameters to pass to the solver.rng
: Random number generator for reproducibility.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.store_measurement
: Whether to store the measurement expectation values. Default isVal(false)
.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-2
andabstol=1e-2
.For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Performance Tip
When sc_ops
contains only a single operator, it is recommended to pass only that operator as the argument. This ensures that the stochastic noise is diagonal, making the simulation faster.
Returns
prob
: TheTimeEvolutionProblem
containing theSDEProblem
for the Stochastic Master Equation time evolution.
QuantumToolbox.smesolveEnsembleProblem Function
smesolveEnsembleProblem(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject,
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
sc_ops::Union{Nothing,AbstractVector,Tuple,AbstractQuantumObject} = nothing;
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
rng::AbstractRNG = default_rng(),
ntraj::Int = 500,
ensemblealg::EnsembleAlgorithm = EnsembleThreads(),
prob_func::Union{Function, Nothing} = nothing,
output_func::Union{Tuple,Nothing} = nothing,
progress_bar::Union{Val,Bool} = Val(true),
store_measurement::Union{Val,Bool} = Val(false),
kwargs...,
)
Generate the SDEProblem for the Stochastic Master Equation time evolution of an open quantum system. This is defined by the following stochastic differential equation:
where
is the Lindblad superoperator, and
Above,
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. It can be either a Ket
or aOperator
.tlist
: List of times at which to save either the state or the expectation values of the system.c_ops
: List of collapse operators. It can be either a Vector
or aTuple
.sc_ops
: List of stochastic collapse operators. It can be either a Vector
, aTuple
or aAbstractQuantumObject
. It is recommended to use the last case when only one operator is provided.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
:NullParameters
of parameters to pass to the solver.rng
: Random number generator for reproducibility.ntraj
: Number of trajectories to use. Default is500
.ensemblealg
: Ensemble method to use. Default toEnsembleThreads()
.prob_func
: Function to use for generating the SDEProblem.output_func
: aTuple
containing theFunction
to use for generating the output of a single trajectory, the (optional)ProgressBar
object, and the (optional)RemoteChannel
object.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.store_measurement
: Whether to store the measurement expectation values. Default isVal(false)
.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-2
andabstol=1e-2
.For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Performance Tip
When sc_ops
contains only a single operator, it is recommended to pass only that operator as the argument. This ensures that the stochastic noise is diagonal, making the simulation faster.
Returns
prob
: TheTimeEvolutionProblem
containing the EnsembleSDEProblem
for the Stochastic Master Equation time evolution.
QuantumToolbox.sesolve Function
sesolve(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject{KetQuantumObject},
tlist::AbstractVector;
alg::OrdinaryDiffEqAlgorithm = Tsit5(),
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
progress_bar::Union{Val,Bool} = Val(true),
inplace::Union{Val,Bool} = Val(true),
kwargs...,
)
Time evolution of a closed quantum system using the Schrödinger equation:
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. tlist
: List of times at which to save either the state or the expectation values of the system.alg
: The algorithm for the ODE solver. The default isTsit5()
.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
: Parameters to pass to the solver. This argument is usually expressed as aNamedTuple
orAbstractVector
of parameters. For more advanced usage, any custom struct can be used.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.inplace
: Whether to use the inplace version of the ODEProblem. The default isVal(true)
. It is recommended to useVal(true)
for better performance, but it is sometimes necessary to useVal(false)
, for example when performing automatic differentiation using Zygote.jl.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-6
andabstol=1e-8
.For more details about
alg
please refer toDifferentialEquations.jl
(ODE Solvers)For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionSol
: The solution of the time evolution. See alsoTimeEvolutionSol
QuantumToolbox.mesolve Function
mesolve(
H::Union{AbstractQuantumObject,Tuple},
ψ0::QuantumObject,
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
alg::OrdinaryDiffEqAlgorithm = Tsit5(),
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
progress_bar::Union{Val,Bool} = Val(true),
inplace::Union{Val,Bool} = Val(true),
kwargs...,
)
Time evolution of an open quantum system using Lindblad master equation:
where
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. It can be either a Ket
or aOperator
.tlist
: List of times at which to save either the state or the expectation values of the system.c_ops
: List of collapse operators. It can be either a Vector
or aTuple
.alg
: The algorithm for the ODE solver. The default value isTsit5()
.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
: Parameters to pass to the solver. This argument is usually expressed as aNamedTuple
orAbstractVector
of parameters. For more advanced usage, any custom struct can be used.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.inplace
: Whether to use the inplace version of the ODEProblem. The default isVal(true)
. It is recommended to useVal(true)
for better performance, but it is sometimes necessary to useVal(false)
, for example when performing automatic differentiation using Zygote.jl.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-6
andabstol=1e-8
.For more details about
alg
please refer toDifferentialEquations.jl
(ODE Solvers)For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionSol
: The solution of the time evolution. See alsoTimeEvolutionSol
QuantumToolbox.mcsolve Function
mcsolve(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject{KetQuantumObject},
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
alg::OrdinaryDiffEqAlgorithm = Tsit5(),
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
rng::AbstractRNG = default_rng(),
ntraj::Int = 500,
ensemblealg::EnsembleAlgorithm = EnsembleThreads(),
jump_callback::TJC = ContinuousLindbladJumpCallback(),
progress_bar::Union{Val,Bool} = Val(true),
prob_func::Union{Function, Nothing} = nothing,
output_func::Union{Tuple,Nothing} = nothing,
normalize_states::Union{Val,Bool} = Val(true),
kwargs...,
)
Time evolution of an open quantum system using quantum trajectories.
Given a system Hamiltonian
with a non-Hermitian effective Hamiltonian:
To the first-order of the wave function in a small time
where
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
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. tlist
: List of times at which to save either the state or the expectation values of the system.c_ops
: List of collapse operators. It can be either a Vector
or aTuple
.alg
: The algorithm to use for the ODE solver. Default toTsit5()
.e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
: Parameters to pass to the solver. This argument is usually expressed as aNamedTuple
orAbstractVector
of parameters. For more advanced usage, any custom struct can be used.rng
: Random number generator for reproducibility.ntraj
: Number of trajectories to use.ensemblealg
: Ensemble algorithm to use. Default toEnsembleThreads()
.jump_callback
: The Jump Callback type: Discrete or Continuous. The default isContinuousLindbladJumpCallback()
, which is more precise.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.prob_func
: Function to use for generating the ODEProblem.output_func
: aTuple
containing theFunction
to use for generating the output of a single trajectory, the (optional)ProgressBar
object, and the (optional)RemoteChannel
object.normalize_states
: Whether to normalize the states. Default toVal(true)
.kwargs
: The keyword arguments for the ODEProblem.
Notes
ensemblealg
can be one ofEnsembleThreads()
,EnsembleSerial()
,EnsembleDistributed()
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-6
andabstol=1e-8
.For more details about
alg
please refer toDifferentialEquations.jl
(ODE Solvers)For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Returns
sol::TimeEvolutionMCSol
: The solution of the time evolution. See alsoTimeEvolutionMCSol
.
QuantumToolbox.ssesolve Function
ssesolve(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject{KetQuantumObject},
tlist::AbstractVector,
sc_ops::Union{Nothing,AbstractVector,Tuple,AbstractQuantumObject} = nothing;
alg::Union{Nothing,StochasticDiffEqAlgorithm} = nothing,
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
rng::AbstractRNG = default_rng(),
ntraj::Int = 500,
ensemblealg::EnsembleAlgorithm = EnsembleThreads(),
prob_func::Union{Function, Nothing} = nothing,
output_func::Union{Tuple,Nothing} = nothing,
progress_bar::Union{Val,Bool} = Val(true),
store_measurement::Union{Val,Bool} = Val(false),
kwargs...,
)
Stochastic Schrödinger equation evolution of a quantum system given the system Hamiltonian
where
and
Above,
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. tlist
: List of times at which to save either the state or the expectation values of the system.sc_ops
: List of stochastic collapse operators. It can be either a Vector
, aTuple
or aAbstractQuantumObject
. It is recommended to use the last case when only one operator is provided.alg
: The algorithm to use for the stochastic differential equation. Default isSRIW1()
ifsc_ops
is anAbstractQuantumObject
(diagonal noise), andSRA2()
otherwise (non-diagonal noise).e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
:NullParameters
of parameters to pass to the solver.rng
: Random number generator for reproducibility.ntraj
: Number of trajectories to use. Default is500
.ensemblealg
: Ensemble method to use. Default toEnsembleThreads()
.prob_func
: Function to use for generating the SDEProblem.output_func
: aTuple
containing theFunction
to use for generating the output of a single trajectory, the (optional)ProgressBar
object, and the (optional)RemoteChannel
object.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.store_measurement
: Whether to store the measurement results. Default isVal(false)
.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-2
andabstol=1e-2
.For more details about
alg
please refer toDifferentialEquations.jl
(SDE Solvers)For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Performance Tip
When sc_ops
contains only a single operator, it is recommended to pass only that operator as the argument. This ensures that the stochastic noise is diagonal, making the simulation faster.
Returns
sol::TimeEvolutionStochasticSol
: The solution of the time evolution. SeeTimeEvolutionStochasticSol
.
QuantumToolbox.smesolve Function
smesolve(
H::Union{AbstractQuantumObject{OperatorQuantumObject},Tuple},
ψ0::QuantumObject,
tlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
sc_ops::Union{Nothing,AbstractVector,Tuple,AbstractQuantumObject} = nothing;
alg::Union{Nothing,StochasticDiffEqAlgorithm} = nothing,
e_ops::Union{Nothing,AbstractVector,Tuple} = nothing,
params = NullParameters(),
rng::AbstractRNG = default_rng(),
ntraj::Int = 500,
ensemblealg::EnsembleAlgorithm = EnsembleThreads(),
prob_func::Union{Function, Nothing} = nothing,
output_func::Union{Tuple,Nothing} = nothing,
progress_bar::Union{Val,Bool} = Val(true),
store_measurement::Union{Val,Bool} = Val(false),
kwargs...,
)
Stochastic Master Equation time evolution of an open quantum system. This is defined by the following stochastic differential equation:
where
is the Lindblad superoperator, and
Above,
Arguments
H
: Hamiltonian of the system. It can be either a QuantumObject
, aQuantumObjectEvolution
, or aTuple
of operator-function pairs.ψ0
: Initial state of the system. It can be either a Ket
or aOperator
.tlist
: List of times at which to save either the state or the expectation values of the system.c_ops
: List of collapse operators. It can be either a Vector
or aTuple
.sc_ops
: List of stochastic collapse operators. It can be either a Vector
, aTuple
or aAbstractQuantumObject
. It is recommended to use the last case when only one operator is provided.alg
: The algorithm to use for the stochastic differential equation. Default isSRIW1()
ifsc_ops
is anAbstractQuantumObject
(diagonal noise), andSRA2()
otherwise (non-diagonal noise).e_ops
: List of operators for which to calculate expectation values. It can be either aVector
or aTuple
.params
:NullParameters
of parameters to pass to the solver.rng
: Random number generator for reproducibility.ntraj
: Number of trajectories to use. Default is500
.ensemblealg
: Ensemble method to use. Default toEnsembleThreads()
.prob_func
: Function to use for generating the SDEProblem.output_func
: aTuple
containing theFunction
to use for generating the output of a single trajectory, the (optional)ProgressBar
object, and the (optional)RemoteChannel
object.progress_bar
: Whether to show the progress bar. Using non-Val
types might lead to type instabilities.store_measurement
: Whether to store the measurement expectation values. Default isVal(false)
.kwargs
: The keyword arguments for the ODEProblem.
Notes
The states will be saved depend on the keyword argument
saveat
inkwargs
.If
e_ops
is empty, the default value ofsaveat=tlist
(saving the states corresponding totlist
), otherwise,saveat=[tlist[end]]
(only save the final state). You can also specifye_ops
andsaveat
separately.The default tolerances in
kwargs
are given asreltol=1e-2
andabstol=1e-2
.For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Performance Tip
When sc_ops
contains only a single operator, it is recommended to pass only that operator as the argument. This ensures that the stochastic noise is diagonal, making the simulation faster.
Returns
sol::TimeEvolutionStochasticSol
: The solution of the time evolution. SeeTimeEvolutionStochasticSol
.
QuantumToolbox.dfd_mesolve Function
dfd_mesolve(H::Function, ψ0::QuantumObject,
tlist::AbstractVector, c_ops::Function, maxdims::AbstractVector,
dfd_params::NamedTuple=NamedTuple();
alg::OrdinaryDiffEqAlgorithm=Tsit5(),
e_ops::Function=(dim_list) -> Vector{Vector{T1}}([]),
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 toDifferentialEquations.jl
(ODE Solvers)For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
QuantumToolbox.liouvillian Function
liouvillian(H::AbstractQuantumObject, c_ops::Union{Nothing,AbstractVector,Tuple}=nothing, Id_cache=I(prod(H.dimensions)))
Construct the Liouvillian SuperOperator
for a system Hamiltonian
where
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
.
QuantumToolbox.liouvillian_generalized Function
liouvillian_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.
Steady State Solvers
QuantumToolbox.steadystate Function
steadystate(
H::QuantumObject{OpType},
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
solver::SteadyStateSolver = SteadyStateDirectSolver(),
kwargs...,
)
Solve the stationary state based on different solvers.
Parameters
H
: The Hamiltonian or the Liouvillian of the system.c_ops
: The list of the collapse operators.solver
: see documentation Solving for Steady-State Solutions for different solvers.kwargs
: The keyword arguments for the solver.
QuantumToolbox.steadystate_fourier Function
steadystate_fourier(
H_0::QuantumObject,
H_p::QuantumObject,
H_m::QuantumObject,
ωd::Number,
c_ops::Union{Nothing,AbstractVector,Tuple} = nothing;
n_max::Integer = 2,
tol::R = 1e-8,
solver::FSolver = SteadyStateLinearSolver(),
kwargs...,
)
Calculates the steady state of a periodically driven system. Here H_0
is the Hamiltonian or the Liouvillian of the undriven system. Considering a monochromatic drive at frequency H_p
and H_m
, where H_p
oscillates as H_m
oscillates as
SteadyStateLinearSolver
: Solves the linear system of equations.SSFloquetEffectiveLiouvillian
: Solves the effective Liouvillian.
For both cases, n_max
is the number of Fourier components to consider, and tol
is the tolerance for the solver.
In the case of SteadyStateLinearSolver
, the full linear system is solved at once:
This is a tridiagonal linear system of the form
where
and
This will allow to simultaneously obtain all the
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 SteadyStateLinearSolver
returns a list of QuantumObject
, containing the density matrices for each Fourier component (SSFloquetEffectiveLiouvillian
returns only
Note
steadystate_floquet
is a synonym of steadystate_fourier
.
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. H_m::QuantumObject
: The Hamiltonian or the Liouvillian of the part of the drive that oscillates as. ω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 = SteadyStateLinearSolver
: The solver to use.kwargs...
: Additional keyword arguments to be passed to the solver.
QuantumToolbox.SteadyStateDirectSolver Type
SteadyStateDirectSolver()
A solver which solves steadystate
by finding the inverse of Liouvillian SuperOperator
using the standard method given in LinearAlgebra
.
QuantumToolbox.SteadyStateEigenSolver Type
SteadyStateEigenSolver()
A solver which solves steadystate
by finding the zero (or lowest) eigenvalue of Liouvillian SuperOperator
using eigsolve
.
QuantumToolbox.SteadyStateLinearSolver Type
SteadyStateLinearSolver(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
.
Arguments
alg::SciMLLinearSolveAlgorithm=KrylovJL_GMRES()
: algorithms given inLinearSolve.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.
QuantumToolbox.SteadyStateODESolver Type
SteadyStateODESolver(
alg = Tsit5(),
ψ0 = nothing,
tmax = Inf,
)
An ordinary differential equation (ODE) solver for solving steadystate
.
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 true
:
or
Arguments
alg::OrdinaryDiffEqAlgorithm=Tsit5()
: The algorithm to solve the ODE.ψ0::Union{Nothing,QuantumObject}=nothing
: The initial state of the system. If not specified, a random pure state will be generated.tmax::Real=Inf
: The final time step for the steady state problem.
For more details about the solvers, please refer to OrdinaryDiffEq.jl
.
QuantumToolbox.SSFloquetEffectiveLiouvillian Type
SSFloquetEffectiveLiouvillian(steadystate_solver = SteadyStateDirectSolver())
A solver which solves steadystate_fourier
by first extracting an effective time-independent Liouvillian and then using the steadystate_solver
to extract the steadystate..
Parameters
steadystate_solver::SteadyStateSolver=SteadyStateDirectSolver()
: The solver to use for the effective Liouvillian.
Note
This solver is only available for steadystate_fourier
.
Dynamical Shifted Fock method
QuantumToolbox.dsf_mesolve Function
dsf_mesolve(H::Function,
ψ0::QuantumObject,
tlist::AbstractVector, c_ops::Function,
op_list::Vector{TOl},
α0_l::Vector{<:Number}=zeros(length(op_list)),
dsf_params::NamedTuple=NamedTuple();
alg::OrdinaryDiffEqAlgorithm=Tsit5(),
e_ops::Function=(op_list,p) -> Vector{TOl}([]),
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 toDifferentialEquations.jl
(ODE Solvers)For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
QuantumToolbox.dsf_mcsolve Function
dsf_mcsolve(H::Function,
ψ0::QuantumObject,
tlist::AbstractVector, c_ops::Function,
op_list::Vector{TOl},
α0_l::Vector{<:Number}=zeros(length(op_list)),
dsf_params::NamedTuple=NamedTuple();
alg::OrdinaryDiffEqAlgorithm=Tsit5(),
e_ops::Function=(op_list,p) -> Vector{TOl}([]),
params::NamedTuple=NamedTuple(),
δα_list::Vector{<:Real}=fill(0.2, length(op_list)),
ntraj::Int=500,
ensemblealg::EnsembleAlgorithm=EnsembleThreads(),
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 toDifferentialEquations.jl
(ODE Solvers)For more details about
kwargs
please refer toDifferentialEquations.jl
(Keyword Arguments)
Low-rank time evolution
QuantumToolbox.TimeEvolutionLRSol Type
struct TimeEvolutionLRSol
A structure storing the results and some information from solving low-rank master equation 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 intimes
.fexpect::Matrix
: The function values at each time point.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.z::Vector{QuantumObject}
: The `z`` matrix of the low-rank algorithm at each time point.B::Vector{QuantumObject}
: TheB
matrix of the low-rank algorithm at each time point.
QuantumToolbox.lr_mesolveProblem Function
lr_mesolveProblem(
H::QuantumObject{OperatorQuantumObject},
z::AbstractArray{T,2},
B::AbstractArray{T,2},
tlist::AbstractVector,
c_ops::Union{AbstractVector,Tuple}=();
e_ops::Union{AbstractVector,Tuple}=(),
f_ops::Union{AbstractVector,Tuple}=(),
opt::NamedTuple = lr_mesolve_options_default,
kwargs...,
)
Formulates the ODEproblem for the low-rank time evolution of the system. The function is called by lr_mesolve
. For more information about the low-rank master equation, see (Gravina and Savona, 2024).
Arguments
H::QuantumObject
: The Hamiltonian of the system.z::AbstractArray
: The initial z matrix of the low-rank algorithm.B::AbstractArray
: The initial B matrix of the low-rank algorithm.tlist::AbstractVector
: The time steps at which the expectation values and function values are calculated.c_ops::Union{AbstractVector,Tuple}
: The list of the collapse operators.e_ops::Union{AbstractVector,Tuple}
: The list of the operators for which the expectation values are calculated.f_ops::Union{AbstractVector,Tuple}
: The list of the functions for which the function values are calculated.opt::NamedTuple
: The options of the low-rank master equation.kwargs
: Additional keyword arguments.
QuantumToolbox.lr_mesolve Function
lr_mesolve(
H::QuantumObject{OperatorQuantumObject},
z::AbstractArray{T,2},
B::AbstractArray{T,2},
tlist::AbstractVector,
c_ops::Union{AbstractVector,Tuple}=();
e_ops::Union{AbstractVector,Tuple}=(),
f_ops::Union{AbstractVector,Tuple}=(),
opt::NamedTuple = lr_mesolve_options_default,
kwargs...,
)
Time evolution of an open quantum system using the low-rank master equation. For more information about the low-rank master equation, see (Gravina and Savona, 2024).
Arguments
H::QuantumObject
: The Hamiltonian of the system.z::AbstractArray
: The initial z matrix of the low-rank algorithm.B::AbstractArray
: The initial B matrix of the low-rank algorithm.tlist::AbstractVector
: The time steps at which the expectation values and function values are calculated.c_ops::Union{AbstractVector,Tuple}
: The list of the collapse operators.e_ops::Union{AbstractVector,Tuple}
: The list of the operators for which the expectation values are calculated.f_ops::Union{AbstractVector,Tuple}
: The list of the functions for which the function values are calculated.opt::NamedTuple
: The options of the low-rank master equation.kwargs
: Additional keyword arguments.
Correlations and Spectrum
QuantumToolbox.correlation_3op_2t Function
correlation_3op_2t(H::AbstractQuantumObject,
ψ0::Union{Nothing,QuantumObject},
tlist::AbstractVector,
τlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple},
A::QuantumObject,
B::QuantumObject,
C::QuantumObject;
kwargs...)
Returns the two-times correlation function of three operators
If the initial state ψ0
is given as nothing
, then the steadystate
will be used as the initial state. Note that this is only implemented if H
is constant (QuantumObject
).
QuantumToolbox.correlation_3op_1t Function
correlation_3op_1t(H::AbstractQuantumObject,
ψ0::Union{Nothing,QuantumObject},
τlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple},
A::QuantumObject,
B::QuantumObject,
C::QuantumObject;
kwargs...)
Returns the two-time correlation function (with only one time coordinate
If the initial state ψ0
is given as nothing
, then the steadystate
will be used as the initial state. Note that this is only implemented if H
is constant (QuantumObject
).
QuantumToolbox.correlation_2op_2t Function
correlation_2op_2t(H::AbstractQuantumObject,
ψ0::Union{Nothing,QuantumObject},
tlist::AbstractVector,
τlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple},
A::QuantumObject,
B::QuantumObject;
reverse::Bool=false,
kwargs...)
Returns the two-times correlation function of two operators
If the initial state ψ0
is given as nothing
, then the steadystate
will be used as the initial state. Note that this is only implemented if H
is constant (QuantumObject
).
When reverse=true
, the correlation function is calculated as
QuantumToolbox.correlation_2op_1t Function
correlation_2op_1t(H::AbstractQuantumObject,
ψ0::Union{Nothing,QuantumObject},
τlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple},
A::QuantumObject,
B::QuantumObject;
reverse::Bool=false,
kwargs...)
Returns the two-time correlation function (with only one time coordinate
If the initial state ψ0
is given as nothing
, then the steadystate
will be used as the initial state. Note that this is only implemented if H
is constant (QuantumObject
).
When reverse=true
, the correlation function is calculated as
QuantumToolbox.spectrum_correlation_fft Function
spectrum_correlation_fft(tlist, corr; inverse=false)
Calculate the power spectrum corresponding to a two-time correlation function using fast Fourier transform (FFT).
Parameters
tlist::AbstractVector
: List of times at which the two-time correlation function is given.corr::AbstractVector
: List of two-time correlations corresponding to the given time point intlist
.inverse::Bool
: Whether to use the inverse Fourier transform or not. Default tofalse
.
Returns
ωlist
: the list of angular frequencies. Slist
: the list of the power spectrum corresponding to the angular frequencies inωlist
.
QuantumToolbox.spectrum Function
spectrum(H::QuantumObject,
ωlist::AbstractVector,
c_ops::Union{Nothing,AbstractVector,Tuple},
A::QuantumObject{OperatorQuantumObject},
B::QuantumObject{OperatorQuantumObject};
solver::SpectrumSolver=ExponentialSeries(),
kwargs...)
Calculate the spectrum of the correlation function
See also the following list for SpectrumSolver
docstrings:
QuantumToolbox.ExponentialSeries Type
ExponentialSeries(; tol = 1e-14, calc_steadystate = false)
A solver which solves spectrum
by finding the eigen decomposition of the Liouvillian SuperOperator
and calculate the exponential series.
QuantumToolbox.PseudoInverse Type
PseudoInverse(; alg::SciMLLinearSolveAlgorithm = KrylovJL_GMRES())
A solver which solves spectrum
by finding the inverse of Liouvillian SuperOperator
using the alg
orithms given in LinearSolve.jl
.
Entropy and Metrics
QuantumToolbox.entropy_vn Function
entropy_vn(ρ::QuantumObject; base::Int=0, tol::Real=1e-15)
Calculates the Von Neumann entropy
Notes
base
specifies the base of the logarithm to use, and when using the default value0
, the natural logarithm is used.tol
describes the absolute tolerance for detecting the zero-valued eigenvalues of the density matrix.
Examples
Pure state:
julia> ψ = fock(2,0)
Quantum Object: type=Ket dims=[2] size=(2,)
2-element Vector{ComplexF64}:
1.0 + 0.0im
0.0 + 0.0im
julia> entropy_vn(ψ, base=2)
-0.0
Mixed state:
julia> ρ = maximally_mixed_dm(2)
Quantum Object: type=Operator dims=[2] size=(2, 2) ishermitian=true
2×2 Diagonal{ComplexF64, Vector{ComplexF64}}:
0.5-0.0im ⋅
⋅ 0.5-0.0im
julia> entropy_vn(ρ, base=2)
1.0
QuantumToolbox.entropy_relative Function
entropy_relative(ρ::QuantumObject, σ::QuantumObject; base::Int=0, tol::Real=1e-15)
Calculates the quantum relative entropy of
Notes
base
specifies the base of the logarithm to use, and when using the default value0
, the natural logarithm is used.tol
describes the absolute tolerance for detecting the zero-valued eigenvalues of the density matrix.
References
- Nielsen and Chuang (2012), section 11.3.1, page 511
QuantumToolbox.entropy_linear Function
entropy_linear(ρ::QuantumObject)
Calculates the quantum linear entropy
QuantumToolbox.entropy_mutual Function
entropy_mutual(ρAB::QuantumObject, selA, selB; kwargs...)
Calculates the quantum mutual information
Here,
Notes
selA
specifies the indices of the sub-systemA
inρAB.dimensions
. See alsoptrace
.selB
specifies the indices of the sub-systemB
inρAB.dimensions
. See alsoptrace
.kwargs
are the keyword arguments for calculating Von Neumann entropy. See alsoentropy_vn
.
QuantumToolbox.entropy_conditional Function
entropy_conditional(ρAB::QuantumObject, selB; kwargs...)
Calculates the conditional quantum entropy with respect to sub-system
Here,
Notes
selB
specifies the indices of the sub-systemB
inρAB.dimensions
. See alsoptrace
.kwargs
are the keyword arguments for calculating Von Neumann entropy. See alsoentropy_vn
.
QuantumToolbox.entanglement Function
entanglement(ρ::QuantumObject, sel; kwargs...)
Calculates the entanglement entropy by doing the partial trace of ρ
, selecting only the dimensions with the indices contained in the sel
vector, and then use the Von Neumann entropy entropy_vn
.
Notes
ρ
can be either aKet
or anOperator
. But should be a pure state.sel
specifies the indices of the remaining sub-system. See alsoptrace
.kwargs
are the keyword arguments for calculating Von Neumann entropy. See alsoentropy_vn
.
QuantumToolbox.concurrence Function
concurrence(ρ::QuantumObject)
Calculate the concurrence for a two-qubit state.
Notes
QuantumToolbox.negativity Function
negativity(ρ::QuantumObject, subsys::Int; logarithmic::Bool=false)
Compute the negativity
Arguments
ρ::QuantumObject
: The density matrix (ρ.type
must beOperatorQuantumObject
).subsys::Int
: an index that indicates which subsystem to compute the negativity for.logarithmic::Bool
: choose whether to calculate logarithmic negativity or not. Default asfalse
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
0.0 + 0.0im
0.0 + 0.0im
0.7071067811865475 + 0.0im
julia> ρ = ket2dm(Ψ)
Quantum Object: type=Operator dims=[2, 2] size=(4, 4) ishermitian=true
4×4 Matrix{ComplexF64}:
0.5+0.0im 0.0+0.0im 0.0+0.0im 0.5+0.0im
0.0+0.0im 0.0+0.0im 0.0+0.0im 0.0+0.0im
0.0+0.0im 0.0+0.0im 0.0+0.0im 0.0+0.0im
0.5+0.0im 0.0+0.0im 0.0+0.0im 0.5+0.0im
julia> round(negativity(ρ, 2), digits=2)
0.5
QuantumToolbox.fidelity Function
fidelity(ρ::QuantumObject, σ::QuantumObject)
Calculate the fidelity of two QuantumObject
:
Here, the definition is from Nielsen and Chuang (2012). It is the square root of the fidelity defined in Jozsa (1994).
QuantumToolbox.tracedist Function
tracedist(ρ::QuantumObject, σ::QuantumObject)
Calculates the trace distance between two QuantumObject
:
QuantumToolbox.hilbert_dist Function
hilbert_dist(ρ::QuantumObject, σ::QuantumObject)
Calculates the Hilbert-Schmidt distance between two QuantumObject
:
Note that ρ
and σ
must be either Ket
or Operator
.
References
QuantumToolbox.hellinger_dist Function
hellinger_dist(ρ::QuantumObject, σ::QuantumObject)
Calculates the Hellinger distance between two QuantumObject
:
Note that ρ
and σ
must be either Ket
or Operator
.
References
QuantumToolbox.bures_dist Function
bures_dist(ρ::QuantumObject, σ::QuantumObject)
Calculate the Bures distance between two QuantumObject
:
Here, the definition of fidelity
QuantumToolbox.bures_angle Function
bures_angle(ρ::QuantumObject, σ::QuantumObject)
Calculate the Bures angle between two QuantumObject
:
Here, the definition of fidelity
Spin Lattice
QuantumToolbox.Lattice Type
Lattice
A Julia constructor for a lattice object. The lattice object is used to define the geometry of the lattice. Nx
and Ny
are the number of sites in the x and y directions, respectively. N
is the total number of sites. lin_idx
is a LinearIndices
object and car_idx
is a CartesianIndices
object, and they are used to efficiently select sites on the lattice.
QuantumToolbox.multisite_operator Function
multisite_operator(dims::Union{AbstractVector, Tuple}, pairs::Pair{<:Integer,<:QuantumObject}...)
A Julia function for generating a multi-site operator dims
and a list of pairs pairs
where the first element of the pair is the site index and the second element is the operator acting on that site.
Arguments
dims::Union{AbstractVector, Tuple}
: A vector of dimensions of the lattice.pairs::Pair{<:Integer,<:QuantumObject}...
: A list of pairs where the first element of the pair is the site index and the second element is the operator acting on that site.
Returns
QuantumObject
: A QuantumObject
representing the multi-site operator.
Example
julia> op = multisite_operator(Val(8), 5=>sigmax(), 7=>sigmaz());
julia> op.dims
8-element SVector{8, Int64} with indices SOneTo(8):
2
2
2
2
2
2
2
2
QuantumToolbox.DissipativeIsing Function
DissipativeIsing(Jx::Real, Jy::Real, Jz::Real, hx::Real, hy::Real, hz::Real, γ::Real, latt::Lattice; boundary_condition::Union{Symbol, Val} = Val(:periodic_bc), order::Integer = 1)
A Julia constructor for a dissipative Ising model. The function returns the Hamiltonian
and the collapse operators
Arguments
Jx::Real
: The coupling constant in the x-direction.Jy::Real
: The coupling constant in the y-direction.Jz::Real
: The coupling constant in the z-direction.hx::Real
: The magnetic field in the x-direction.hy::Real
: The magnetic field in the y-direction.hz::Real
: The magnetic field in the z-direction.γ::Real
: The local dissipation rate.latt::Lattice
: ALattice
object that defines the geometry of the lattice.boundary_condition::Union{Symbol, Val}
: The boundary conditions of the lattice. The possible inputs areperiodic_bc
andopen_bc
, for periodic or open boundary conditions, respectively. The default value isVal(:periodic_bc)
.order::Integer
: The order of the nearest-neighbour sites. The default value is 1.
Symmetries and Block Diagonalization
QuantumToolbox.block_diagonal_form Function
block_diagonal_form(A::QuantumObject)
Return the block-diagonal form of a QuantumObject
. This is very useful in the presence of symmetries.
Arguments
A::QuantumObject
: The quantum object.
Returns
The BlockDiagonalForm
of A
.
QuantumToolbox.BlockDiagonalForm Type
struct BlockDiagonalForm
A type for storing a block-diagonal form of a matrix.
Fields
B::DT
: The block-diagonal matrix. It can be a sparse matrix or aQuantumObject
.P::DT
: The permutation matrix. It can be a sparse matrix or aQuantumObject
.blocks::AbstractVector
: The blocks of the block-diagonal matrix.block_sizes::AbstractVector
: The sizes of the blocks.
Miscellaneous
QuantumToolbox.wigner Function
wigner(
state::QuantumObject{OpType},
xvec::AbstractVector,
yvec::AbstractVector;
g::Real = √2,
method::WignerSolver = WignerClenshaw(),
)
Generates the Wigner quasipropability distribution of state
at points xvec + 1im * yvec
in phase space. The g
parameter is a scaling factor related to the value of
The method
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
method has an optional parallel
parameter which defaults to true
and uses multithreading to speed up the calculation.
Arguments
state::QuantumObject
: The quantum state for which the Wigner function is calculated. It can be either aKetQuantumObject
,BraQuantumObject
, orOperatorQuantumObject
.xvec::AbstractVector
: The x-coordinates of the phase space grid.yvec::AbstractVector
: The y-coordinates of the phase space grid.g::Real
: The scaling factor related to the value ofin the commutation relation via . method::WignerSolver
: The method used to calculate the Wigner function. It can be eitherWignerLaguerre()
orWignerClenshaw()
, withWignerClenshaw()
as default. TheWignerLaguerre
method has the optionalparallel
andtol
parameters, with default valuestrue
and1e-14
, respectively.
Returns
W::Matrix
: The Wigner function of the state at the pointsxvec + 1im * yvec
in phase space.
Example
julia> ψ = fock(10, 0) + fock(10, 1) |> normalize
Quantum Object: type=Ket dims=[10] size=(10,)
10-element Vector{ComplexF64}:
0.7071067811865475 + 0.0im
0.7071067811865475 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
0.0 + 0.0im
julia> xvec = range(-5, 5, 200)
-5.0:0.05025125628140704:5.0
julia> wig = wigner(ψ, xvec, xvec);
or taking advantage of the parallel computation of the WignerLaguerre
method
julia> ρ = ket2dm(ψ) |> to_sparse;
julia> wig = wigner(ρ, xvec, xvec, method=WignerLaguerre(parallel=true));
Linear Maps
QuantumToolbox.AbstractLinearMap Type
AbstractLinearMap{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 abstract type helps to define this map when the matrix is not explicitly available.
Methods
Base.eltype(A)
: Returns the element typeT
.Base.size(A)
: Returns the sizeA.size
.Base.size(A, i)
: Returns thei
-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
end
function LinearAlgebra.mul!(y::AbstractVector, A::ArnoldiLindbladIntegratorMap, x::AbstractVector)
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.
Utility functions
QuantumToolbox.versioninfo Function
QuantumToolbox.versioninfo(io::IO=stdout)
Command line output of information on QuantumToolbox, dependencies, and system information, same as QuantumToolbox.about
.
QuantumToolbox.about Function
QuantumToolbox.about(io::IO=stdout)
Command line output of information on QuantumToolbox, dependencies, and system information, same as QuantumToolbox.versioninfo
.
QuantumToolbox.gaussian Function
gaussian(x::Number, μ::Number, σ::Number)
Returns the gaussian function
QuantumToolbox.n_thermal Function
n_thermal(ω::Real, ω_th::Real)
Return the number of photons in thermal equilibrium for an harmonic oscillator mode with frequency
where
QuantumToolbox.PhysicalConstants Constant
const 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 unitG
: Newtonian constant of gravitation with unith
: (exact) Planck constant with unitħ
: reduced Planck constant with unite
: (exact) elementary charge with unitμ0
: vacuum magnetic permeability with unitϵ0
: vacuum electric permittivity with unitk
: (exact) Boltzmann constant with unitNA
: (exact) Avogadro constant with unit
Examples
julia> PhysicalConstants.ħ
1.0545718176461565e-34
QuantumToolbox.convert_unit Function
convert_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:GHz
: Giga-Hertz multiplied by Planck constant:K
: Kelvin multiplied by Boltzmann constant:mK
: milli-Kelvin multiplied by Boltzmann constant
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> round(convert_unit(1, :meV, :mK), digits=4)
11604.5181
QuantumToolbox.row_major_reshape Function
row_major_reshape(Q::AbstractArray, shapes...)
Reshapes Q
in the row-major order, as numpy.
QuantumToolbox.meshgrid Function
meshgrid(x::AbstractVector, y::AbstractVector)
Equivalent to numpy meshgrid.
Visualization
QuantumToolbox.plot_wigner Function
plot_wigner(
state::QuantumObject{OpType};
library::Union{Val,Symbol}=Val(:CairoMakie),
kwargs...
) where {OpType<:Union{BraQuantumObject,KetQuantumObject,OperatorQuantumObject}
Plot the Wigner quasipropability distribution of state
using the wigner
function.
The library
keyword argument specifies the plotting library to use, defaulting to CairoMakie
.
Arguments
state::QuantumObject
: The quantum state for which to plot the Wigner distribution.library::Union{Val,Symbol}
: The plotting library to use. Default isVal(:CairoMakie)
.kwargs...
: Additional keyword arguments to pass to the plotting function. See the documentation for the specific plotting library for more information.
Import library first
The plotting libraries must first be imported before using them with this function.
Beware of type-stability!
If you want to keep type stability, it is recommended to use Val(:CairoMakie)
instead of :CairoMakie
as the plotting library. See this link and the related Section about type stability for more details.
plot_wigner(
library::Val{:CairoMakie},
state::QuantumObject{OpType};
xvec::Union{Nothing,AbstractVector} = nothing,
yvec::Union{Nothing,AbstractVector} = nothing,
g::Real = √2,
method::WignerSolver = WignerClenshaw(),
projection::Union{Val,Symbol} = Val(:two_dim),
location::Union{GridPosition,Nothing} = nothing,
colorbar::Bool = false,
kwargs...
) where {OpType}
Plot the Wigner quasipropability distribution of state
using the CairoMakie
plotting library.
Arguments
library::Val{:CairoMakie}
: The plotting library to use.state::QuantumObject
: The quantum state for which the Wigner function is calculated. It can be either aKet
,Bra
, orOperator
.xvec::AbstractVector
: The x-coordinates of the phase space grid. Defaults to a linear range from -7.5 to 7.5 with 200 points.yvec::AbstractVector
: The y-coordinates of the phase space grid. Defaults to a linear range from -7.5 to 7.5 with 200 points.g::Real
: The scaling factor related to the value ofin the commutation relation via . method::WignerSolver
: The method used to calculate the Wigner function. It can be eitherWignerLaguerre()
orWignerClenshaw()
, withWignerClenshaw()
as default. TheWignerLaguerre
method has the optionalparallel
andtol
parameters, with default valuestrue
and1e-14
, respectively.projection::Union{Val,Symbol}
: Whether to plot the Wigner function in 2D or 3D. It can be eitherVal(:two_dim)
orVal(:three_dim)
, withVal(:two_dim)
as default.location::Union{GridPosition,Nothing}
: The location of the plot in the layout. Ifnothing
, the plot is created in a new figure. Default isnothing
.colorbar::Bool
: Whether to include a colorbar in the plot. Default isfalse
.kwargs...
: Additional keyword arguments to pass to the plotting function.
Returns
fig
: The figure object.ax
: The axis object.hm
: Either the heatmap or surface object, depending on the projection.
Import library first
CairoMakie
must first be imported before using this function.
Beware of type-stability!
If you want to keep type stability, it is recommended to use Val(:two_dim)
and Val(:three_dim)
instead of :two_dim
and :three_dim
, respectively. Also, specify the library as Val(:CairoMakie)
See this link and the related Section about type stability for more details.