API
Hilbert Spaces
QuantumCumulants.HilbertSpace
— TypeHilbertSpace
Abstract type for representing Hilbert spaces.
QuantumCumulants.ProductSpace
— TypeProductSpace <: HilbertSpace
Stores a composite HilbertSpace
consisting of multiple subspaces. Generally created by computing the tensor product ⊗
of subspaces.
QuantumCumulants.FockSpace
— TypeFockSpace <: HilbertSpace
HilbertSpace
defining a Fock space for bosonic operators. See also: Destroy
, Create
QuantumCumulants.NLevelSpace
— TypeNLevelSpace <: HilbertSpace
NLevelSpace(name::Symbol,levels,GS=levels[1])
Define a HilbertSpace
for an object consisting of N
discrete energy levels. The given levels
must be an integer specifying the number of levels, or an iterable collection of levels. The argument GS
specifies which state should be treated as ground state and is rewritten using population conservation during simplification. See also: Transition
Examples:
julia> ha = NLevelSpace(:a,3)
ℋ(a)
julia> ha = NLevelSpace(:a,(:g,:e))
ℋ(a)
QuantumCumulants.ClusterSpace
— TypeClusterSpace <: HilbertSpace
ClusterSpace(original_space,N,order)
A Hilbert space representing N
identical copies of another Hilbert space, with correlations up to a specified order
.
QuantumCumulants.:⊗
— Function⊗(spaces::HilbertSpace...)
Create a ProductSpace
consisting of multiple subspaces. Unicode \otimes<tab>
alias of tensor
Examples
julia> hf = FockSpace(:f)
ℋ(f)
julia> ha = NLevelSpace(:a,2)
ℋ(a)
julia> h = hf⊗ha
ℋ(f) ⊗ ℋ(a)
QuantumCumulants.tensor
— Functiontensor(spaces::HilbertSpace...)
Create a ProductSpace
consisting of multiple subspaces. See also ⊗
.
q-Numbers
QuantumCumulants.QSym
— TypeQSym <: QNumber
Abstract type representing fundamental operator types.
QuantumCumulants.QTerm
— TypeQTerm <: QNumber
Abstract type representing noncommutative expressions.
QuantumCumulants.@qnumbers
— Macro@qnumbers
Convenience macro for the construction of operators.
Examples
julia> h = FockSpace(:fock)
ℋ(fock)
julia> @qnumbers a::Destroy(h)
(a,)
julia> h = FockSpace(:one) ⊗ FockSpace(:two)
ℋ(one) ⊗ ℋ(two)
julia> @qnumbers b::Destroy(h,2)
(b,)
QuantumCumulants.Destroy
— TypeDestroy <: QSym
Bosonic operator on a FockSpace
representing the quantum harmonic oscillator annihilation operator.
QuantumCumulants.Create
— TypeCreate <: QSym
Bosonic operator on a FockSpace
representing the quantum harmonic oscillator creation operator.
QuantumCumulants.Transition
— TypeTransition <: QSym
Transition(h::NLevelSpace,name::Symbol,i,j)
Fundamental operator defining a transition from level j
to level i
on a NLevelSpace
. The notation corresponds to Dirac notation, i.e. the above is equivalent to |i⟩⟨j|
.
Examples
julia> ha = NLevelSpace(:a,(:g,:e))
ℋ(a)
julia> σ = Transition(ha,:σ,:g,:e)
σge
Mean field
QuantumCumulants.meanfield
— Functionmeanfield(ops::Vector,H::QNumber)
meanfield(op::QNumber,H::QNumber)
meanfield(ops::Vector,H::QNumber,J::Vector;
Jdagger::Vector=adjoint.(J),rates=ones(length(J)))
meanfield(op::QNumber,H::QNumber,J::Vector;
Jdagger::Vector=adjoint.(J),rates=ones(length(J)))
Compute the set of equations for the operators in ops
under the Hamiltonian H
and with loss operators contained in J
. The resulting equation is equivalent to the Quantum-Langevin equation where noise is neglected.
Arguments
*ops::Vector
: The operators of which the equations are to be computed. *H::QNumber
: The Hamiltonian describing the reversible dynamics of the system. *J::Vector{<:QNumber}
: A vector containing the collapse operators of the system. A term of the form $\sum_i J_i^\dagger O J_i - \frac{1}{2}\left(J_i^\dagger J_i O + OJ_i^\dagger J_i\right)$ is added to the Heisenberg equation.
Optional argumentes
*Jdagger::Vector=adjoint.(J)
: Vector containing the hermitian conjugates of the collapse operators. *rates=ones(length(J))
: Decay rates corresponding to the collapse operators in J
. *multithread=false
: Specify whether the derivation of equations for all operators in ops
should be multithreaded using Threads.@threads
. *simplify=true
: Specify whether the derived equations should be simplified. *order=nothing
: Specify to which order
a cumulant_expansion
is performed. If nothing
, this step is skipped. *mix_choice=maximum
: If the provided order
is a Vector
, mix_choice
determines which order
to prefer on terms that act on multiple Hilbert spaces. *iv=SymbolicUtils.Sym{Real}(:t)
: The independent variable (time parameter) of the system.
QuantumCumulants.commutator
— Functioncommutator(a,b)
Computes the commutator a*b - b*a
.
QuantumCumulants.acts_on
— Functionacts_on(op)
Shows on which Hilbert space op
acts. For QSym
types, this returns an Integer, whereas for a Term
it returns a Vector{Int}
whose entries specify all subspaces on which the expression acts.
QuantumCumulants.MeanfieldEquations
— TypeMeanfieldEquations <: AbstractMeanfieldEquations
Type defining a system of differential equations, where lhs
is a vector of derivatives and rhs
is a vector of expressions. In addition, it keeps track of the Hamiltonian, the collapse operators and the corresponding decay rates of the system.
Fields
*equations
: Vector of the differential equations of averages. *operator_equations
: Vector of the operator differential equations. *states
: Vector containing the averages on the left-hand-side of the equations. *operators
: Vector containing the operators on the left-hand-side of the equations. *hamiltonian
: Operator defining the system Hamiltonian. *jumps
: Vector of operators specifying the decay processes. *jumps
: Vector of operators specifying the adjoint of the decay processes. *rates
: Decay rates corresponding to the jumps
. *iv
: The independent variable (time parameter) of the system. *varmap
: Vector of pairs that map the averages to time-dependent variables. That format is necessary for ModelingToolkit functionality. *order
: The order at which the cumulant_expansion
has been performed.
Symbolic Numbers
QuantumCumulants.CNumber
— TypeQuantumCumulants.Parameter
— TypeParameter <: CNumber
Type used as symbolic type in a SymbolicUtils.Sym
variable to represent a parameter.
QuantumCumulants.cnumbers
— Functioncnumbers(symbols::Symbol...)
cnumbers(s::String)
Create symbolic cnumbers.
Expamples
julia> ps = cnumbers(:a, :b)
(a, b)
julia> cnumbers("a b") == ps
true
QuantumCumulants.@cnumbers
— Macro@cnumbers(ps...)
Convenience macro to quickly define symbolic cnumbers.
Examples
julia> @cnumbers ω κ
(ω, κ)
Average
QuantumCumulants.average
— Functionaverage(::QNumber)
average(::QNumber,order)
Compute the average of an operator. If order
is given, the cumulant_expansion
up to that order is computed immediately.
QuantumCumulants.cumulant_expansion
— Functioncumulant_expansion(avg, order::Int)
For an average
of an operator, expand it in terms of moments up to order
neglecting their joint cumulant.
See also: https://en.wikipedia.org/wiki/Cumulant#Joint_cumulants
Examples
julia> avg = average(a*b)
⟨a*b⟩
julia> cumulant_expansion(avg,1)
(⟨a⟩*⟨b⟩)
julia> avg = average(a*b*c)
⟨a*b*c⟩
julia> cumulant_expansion(avg,2)
((⟨a*b⟩*⟨c⟩)+(⟨a*c⟩*⟨b⟩)+(⟨a⟩*⟨b*c⟩)+(-2*⟨a⟩*⟨b⟩*⟨c⟩))
Optional arguments
*simplify=true: Specify whether the result should be simplified. *kwargs...: Further keyword arguments being passed to simplification.
QuantumCumulants.cumulant
— Functioncumulant(x,n=get_order(x);simplify=true,kwargs...)
Compute the n
th cumulant of x
(either an operator or an average). The output is simplified when simplify=true
. Further keyword arguments are passed on to simplification.
Examples
julia> cumulant(a)
⟨a⟩
julia> cumulant(a*b)
(⟨a*b⟩+(-1*⟨a⟩*⟨b⟩))
julia> cumulant(a*b,1)
⟨a*b⟩
julia> cumulant(a*b,3)
0
QuantumCumulants.get_order
— Functionget_order(arg)
Compute the order of a given argument. This is the order used to decide whether something should be expanded using a cumulant_expansion
method.
Examples
julia> get_order(a)
1
julia> get_order(a*b)
2
julia> get_order(1)
0
Correlation functions
QuantumCumulants.CorrelationFunction
— Typestruct CorrelationFunction
Type representing the two-time first-order correlation function of two operators.
QuantumCumulants.Spectrum
— Typestruct Spectrum
Type representing the spectrum, i.e. the Fourier transform of a CorrelationFunction
in steady state.
To actually compute the spectrum at a frequency ω
, construct the type on top of a correlation function and call it with Spectrum(c)(ω,usteady,p0)
.
QuantumCumulants.correlation_u0
— Functioncorrelation_u0(c::CorrelationFunction, u_end)
Find the vector containing the correct initial values when numerical solving the time evolution for the correlation function.
See also: CorrelationFunction
correlation_p0
QuantumCumulants.correlation_p0
— Functioncorrelation_p0(c::CorrelationFunction, u_end, ps=Pair[])
Find all occurring steady-state values and add them to a list of parameters to pass this to the ODEProblem
.
See also: CorrelationFunction
correlation_u0
Utility functions
QuantumCumulants.find_missing
— Functionfind_missing(me::MeanfieldEquations, vs_adj=nothing, get_adjoints=true)
Find all averages on the right-hand-side of in me.equations
that are not listed me.states
. For a complete system this list is empty.
Optional arguments
*vs_adj
: List of the complex conjugates of me.states
. If set to nothing
the list is generated internally. *get_adjoints=true
: Specify whether a complex conjugate of an average should be explicitly listed as missing.
QuantumCumulants.find_operators
— Functionfind_operators(::HilbertSpace, order; names=nothing)
Find all operators that fully define a system up to the given order
.
QuantumCumulants.complete
— Functioncomplete(de::MeanfieldEquations)
From a set of differential equation of averages, find all averages that are missing and derive the corresponding equations of motion. Uses find_missing
and meanfield
to do so.
Optional arguments
*order=de.order
: The order at which the cumulant_expansion
is performed on the newly derived equations. If nothing
, the order is inferred from the existing equations. *filter_func=nothing
: Custom function that specifies whether some averages should be ignored when completing a system. This works by calling filter!(filter_func, missed)
where missed
is the vector resulting from find_missing
. Occurrences of averages for which filter_func
returns false
are substituted to 0. *kwargs...
: Further keyword arguments are passed on to meanfield
and simplification.
see also: find_missing
, meanfield
QuantumCumulants.complete!
— Functioncomplete!(de::MeanfieldEquations)
In-place version of complete
QuantumCumulants.unique_ops
— Functionunique_ops(ops)
For a given list of operators, return only unique ones taking into account their adjoints.
QuantumCumulants.unique_ops!
— Functionunique_ops!(ops)
In-place version of unique_ops
.
QuantumCumulants.fundamental_operators
— Functionfundamental_operators(::HilbertSpace)
Return all fundamental operators for a given Hilbertspace. For example, a FockSpace
only has one fundamental operator, Destroy
.
QuantumCumulants.transition_superscript
— Functiontransition_superscript(::Bool)
Specify whether the indices in a Transition
operator should be printed as superscript. Default is true
. If set to false
, the indices corresponding to the levels are printed as subscript.