API
Quantum object functions
QuPhys.BraQuantumObject
— TypeBraQuantumObject <: QuantumObjectType
Abstract type representing a bra state $\bra{\psi}$.
QuPhys.KetQuantumObject
— TypeKetQuantumObject <: QuantumObjectType
Abstract type representing a ket state $\ket{\psi}$.
QuPhys.OperatorQuantumObject
— TypeOperatorQuantumObject <: QuantumObjectType
Abstract type representing an operator $\hat{O}$.
QuPhys.SuperOperatorQuantumObject
— TypeSuperOperatorQuantumObject <: QuantumObjectType
Abstract type representing a super-operator $\hat{\mathcal{O}}$.
QuPhys.QuantumObject
— Typemutable struct QuantumObject{MT<:AbstractArray,ObjType<:QuantumObjectType}
+API · QuPhys.jl API
Quantum object functions
QuPhys.BraQuantumObject
— TypeBraQuantumObject <: QuantumObjectType
Abstract type representing a bra state $\bra{\psi}$.
sourceQuPhys.KetQuantumObject
— TypeKetQuantumObject <: QuantumObjectType
Abstract type representing a ket state $\ket{\psi}$.
sourceQuPhys.OperatorQuantumObject
— TypeOperatorQuantumObject <: QuantumObjectType
Abstract type representing an operator $\hat{O}$.
sourceQuPhys.SuperOperatorQuantumObject
— TypeSuperOperatorQuantumObject <: QuantumObjectType
Abstract type representing a super-operator $\hat{\mathcal{O}}$.
sourceQuPhys.QuantumObject
— Typemutable struct QuantumObject{MT<:AbstractArray,ObjType<:QuantumObjectType}
data::MT
type::Type{ObjType}
dims::Vector{Int}
@@ -13,7 +13,7 @@
⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢
julia> a isa QuantumObject
-true
sourceQuPhys.ket2dm
— Functionket2dm(ψ::QuantumObject)
Transform the ket state $\ket{\psi}$ into a pure density matrix $\hat{\rho} = \dyad{\psi}$.
sourceQuPhys.isbra
— Functionisbra(A::QuantumObject)
Checks if the QuantumObject
A
is a BraQuantumObject
state.
sourceQuPhys.isket
— Functionisket(A::QuantumObject)
Checks if the QuantumObject
A
is a KetQuantumObject
state.
sourceQuPhys.isoper
— Functionisoper(A::QuantumObject)
Checks if the QuantumObject
A
is a OperatorQuantumObject
state.
sourceQuPhys.issuper
— Functionissuper(A::QuantumObject)
Checks if the QuantumObject
A
is a SuperOperatorQuantumObject
state.
sourceBase.size
— Functionsize(A::QuantumObject)
Returns the size of the matrix or vector corresponding to the QuantumObject
A
.
sourceBase.length
— Functionlength(A::QuantumObject)
Returns the length of the matrix or vector corresponding to the QuantumObject
A
.
sourceLinearAlgebra.tr
— Functiontr(A::QuantumObject})
Returns the trace of A
.
Examples
julia> a = destroy(20)
+true
sourceQuPhys.ket2dm
— Functionket2dm(ψ::QuantumObject)
Transform the ket state $\ket{\psi}$ into a pure density matrix $\hat{\rho} = \dyad{\psi}$.
sourceQuPhys.isbra
— Functionisbra(A::QuantumObject)
Checks if the QuantumObject
A
is a BraQuantumObject
state.
sourceQuPhys.isket
— Functionisket(A::QuantumObject)
Checks if the QuantumObject
A
is a KetQuantumObject
state.
sourceQuPhys.isoper
— Functionisoper(A::QuantumObject)
Checks if the QuantumObject
A
is a OperatorQuantumObject
state.
sourceQuPhys.issuper
— Functionissuper(A::QuantumObject)
Checks if the QuantumObject
A
is a SuperOperatorQuantumObject
state.
sourceBase.size
— Functionsize(A::QuantumObject)
Returns the size of the matrix or vector corresponding to the QuantumObject
A
.
sourceBase.length
— Functionlength(A::QuantumObject)
Returns the length of the matrix or vector corresponding to the QuantumObject
A
.
sourceLinearAlgebra.tr
— Functiontr(A::QuantumObject})
Returns the trace of A
.
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:
⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀
@@ -23,7 +23,7 @@
⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢
julia> tr(a' * a)
-190.0 + 0.0im
sourceLinearAlgebra.norm
— Functionnorm(A::QuantumObject)
Returns the norm of A
.
Examples
julia> ψ = fock(10, 2)
+190.0 + 0.0im
sourceLinearAlgebra.norm
— Functionnorm(A::QuantumObject)
Returns the norm of A
.
Examples
julia> ψ = fock(10, 2)
Quantum Object: type=Ket dims=[10] size=(10,)
10-element Vector{ComplexF64}:
0.0 + 0.0im
@@ -38,7 +38,7 @@
0.0 + 0.0im
julia> norm(ψ)
-1.0
sourceBase.kron
— Functionkron(A::QuantumObject, B::QuantumObject)
Returns the Kronecker product $\hat{A} \otimes \hat{B}$.
Examples
julia> a = destroy(20)
+1.0
sourceBase.kron
— Functionkron(A::QuantumObject, B::QuantumObject)
Returns the Kronecker product $\hat{A} \otimes \hat{B}$.
Examples
julia> a = destroy(20)
Quantum Object: type=Operator dims=[20] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 19 stored entries:
⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀
@@ -70,7 +70,7 @@
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡀⠀⠀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡀⠀⠀⠀
⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⢦⡀⠀
-⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⠦
sourceGeneral functions
QuPhys.row_major_reshape
— Functionrow_major_reshape(Q::AbstractArray, shapes...)
Reshapes Q
in the row-major order, as numpy.
sourceQuPhys.meshgrid
— Functionmeshgrid(x::AbstractVector, y::AbstractVector)
Equivalent to numpy meshgrid.
sourceQuPhys.sparse_to_dense
— Functionsparse_to_dense(A::QuantumObject)
Converts a sparse QuantumObject to a dense QuantumObject.
sourceQuPhys.dense_to_sparse
— Functiondense_to_sparse(A::QuantumObject)
Converts a dense QuantumObject to a sparse QuantumObject.
sourceQuPhys.tidyup
— Functiontidyup(A::QuantumObject, tol::Real=1e-14)
Removes those elements of a QuantumObject A
whose absolute value is less than tol
.
sourceQuPhys.tidyup!
— Functiontidyup!(A::QuantumObject, tol::Real=1e-14)
Removes those elements of a QuantumObject A
whose absolute value is less than tol
. In-place version of tidyup
.
sourceQuPhys.gaussian
— Functiongaussian(x::Number, μ::Number, σ::Number)
Returns the gaussian function $\exp \left[- \frac{(x - \mu)^2}{2 \sigma^2} \right]$, where $\mu$ and $\sigma^2$ are the mean and the variance respectively.
sourceQuPhys.ptrace
— Functionptrace(QO::QuantumObject, sel::Vector{Int})
Partial trace of a quantum state QO
leaving only the dimensions with the indices present in the sel
vector.
Examples
Two qubits in the state $\ket{\psi} = \ket{e,g}$:
julia> ψ = kron(fock(2,0), fock(2,1))
+⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠀⠙⠦
sourceGeneral functions
QuPhys.row_major_reshape
— Functionrow_major_reshape(Q::AbstractArray, shapes...)
Reshapes Q
in the row-major order, as numpy.
sourceQuPhys.meshgrid
— Functionmeshgrid(x::AbstractVector, y::AbstractVector)
Equivalent to numpy meshgrid.
sourceQuPhys.sparse_to_dense
— Functionsparse_to_dense(A::QuantumObject)
Converts a sparse QuantumObject to a dense QuantumObject.
sourceQuPhys.dense_to_sparse
— Functiondense_to_sparse(A::QuantumObject)
Converts a dense QuantumObject to a sparse QuantumObject.
sourceQuPhys.tidyup
— Functiontidyup(A::QuantumObject, tol::Real=1e-14)
Removes those elements of a QuantumObject A
whose absolute value is less than tol
.
sourceQuPhys.tidyup!
— Functiontidyup!(A::QuantumObject, tol::Real=1e-14)
Removes those elements of a QuantumObject A
whose absolute value is less than tol
. In-place version of tidyup
.
sourceQuPhys.gaussian
— Functiongaussian(x::Number, μ::Number, σ::Number)
Returns the gaussian function $\exp \left[- \frac{(x - \mu)^2}{2 \sigma^2} \right]$, where $\mu$ and $\sigma^2$ are the mean and the variance respectively.
sourceQuPhys.ptrace
— Functionptrace(QO::QuantumObject, sel::Vector{Int})
Partial trace of a quantum state QO
leaving only the dimensions with the indices present in the sel
vector.
Examples
Two qubits in the state $\ket{\psi} = \ket{e,g}$:
julia> ψ = kron(fock(2,0), fock(2,1))
Quantum Object: type=Ket dims=[2, 2] size=(4,)
4-element Vector{ComplexF64}:
0.0 + 0.0im
@@ -94,7 +94,7 @@
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
sourceQuPhys.entropy_vn
— Functionentropy_vn(ρ::QuantumObject; base::Int=0, tol::Real=1e-15)
Calculates the Von Neumann entropy $S = - \Tr \left[ \hat{\rho} \log \left( \hat{\rho} \right) \right]$ where $\hat{\rho}$ is the density matrix of the system.
The base
parameter specifies the base of the logarithm to use, and when using the default value 0, the natural logarithm is used. The tol
parameter describes the absolute tolerance for detecting the zero-valued eigenvalues of the density matrix $\hat{\rho}$.
Examples
Pure state:
julia> ψ = fock(2,0)
+ 0.0+0.0im 0.5+0.0im
sourceQuPhys.entropy_vn
— Functionentropy_vn(ρ::QuantumObject; base::Int=0, tol::Real=1e-15)
Calculates the Von Neumann entropy $S = - \Tr \left[ \hat{\rho} \log \left( \hat{\rho} \right) \right]$ where $\hat{\rho}$ is the density matrix of the system.
The base
parameter specifies the base of the logarithm to use, and when using the default value 0, the natural logarithm is used. The tol
parameter describes the absolute tolerance for detecting the zero-valued eigenvalues of the density matrix $\hat{\rho}$.
Examples
Pure state:
julia> ψ = fock(2,0)
Quantum Object: type=Ket dims=[2] size=(2,)
2-element Vector{ComplexF64}:
1.0 + 0.0im
@@ -114,13 +114,13 @@
0.0+0.0im 0.5+0.0im
julia> entropy_vn(ρ, base=2)
-1.0
sourceQuPhys.entanglement
— Functionentanglement(QO::QuantumObject, sel::Vector)
Calculates the entanglement by doing the partial trace of QO
, selecting only the dimensions with the indices contained in the sel
vector, and then using the Von Neumann entropy entropy_vn
.
sourceQuPhys.expect
— Functionexpect(O::QuantumObject, ψ::QuantumObject)
Expectation value of the operator O
with the state ψ
. The latter can be a KetQuantumObject
, BraQuantumObject
or a OperatorQuantumObject
. If ψ
is a density matrix, the function calculates $\Tr \left[ \hat{O} \hat{\psi} \right]$, while if ψ
is a state, the function calculates $\mel{\psi}{\hat{O}}{\psi}$.
The function returns a real number if the operator is hermitian, and returns a complex number otherwise.
Examples
julia> ψ = 1 / √2 * (fock(10,2) + fock(10,4));
+1.0
sourceQuPhys.entanglement
— Functionentanglement(QO::QuantumObject, sel::Vector)
Calculates the entanglement by doing the partial trace of QO
, selecting only the dimensions with the indices contained in the sel
vector, and then using the Von Neumann entropy entropy_vn
.
sourceQuPhys.expect
— Functionexpect(O::QuantumObject, ψ::QuantumObject)
Expectation value of the operator O
with the state ψ
. The latter can be a KetQuantumObject
, BraQuantumObject
or a OperatorQuantumObject
. If ψ
is a density matrix, the function calculates $\Tr \left[ \hat{O} \hat{\psi} \right]$, while if ψ
is a state, the function calculates $\mel{\psi}{\hat{O}}{\psi}$.
The function returns a real number if the operator is hermitian, and returns a complex number otherwise.
Examples
julia> ψ = 1 / √2 * (fock(10,2) + fock(10,4));
julia> a = destroy(10);
julia> expect(a' * a, ψ) ≈ 3
-true
sourceQuPhys.wigner
— Functionwigner(state::QuantumObject, xvec::AbstractVector, yvec::AbstractVector; g::Real=√2,
- solver::WignerSolver=WignerLaguerre())
Generates the Wigner quasipropability distribution of state
at points xvec + 1im * yvec
. The g
parameter is a scaling factor related to the value of $\hbar$ in the commutation relation $[x, y] = i \hbar$ via $\hbar=2/g^2$ giving the default value $\hbar=1$.
The solver
parameter can be either WignerLaguerre()
or WignerClenshaw()
. The former uses the Laguerre polynomial expansion of the Wigner function, while the latter uses the Clenshaw algorithm. The Laguerre expansion is faster for sparse matrices, while the Clenshaw algorithm is faster for dense matrices. The WignerLaguerre
solver has an optional parallel
parameter which defaults to true
and uses multithreading to speed up the calculation.
sourceQuPhys.get_coherence
— Functionget_coherence(ψ::QuantumObject)
Get the coherence value $\alpha$ by measuring the expectation value of the destruction operator $\hat{a}$ on the state $\ket{\psi}$.
It returns both $\alpha$ and the state $\ket{\delta_\psi} = \exp ( \bar{\alpha} \hat{a} - \alpha \hat{a}^\dagger )$. The latter corresponds to the quantum fulctuations around the coherent state $\ket{\alpha}$.
sourceQuPhys.n_th
— Functionn_th(ω::Number, T::Real)
Gives the mean number of excitations in a mode with frequency ω at temperature T: $n_{\rm th} (\omega, T) = \frac{1}{e^{\omega/T} - 1}$
sourceQuPhys.get_data
— Functionget_data(A::QuantumObject)
Returns the data of a QuantumObject.
sourceQuPhys.mat2vec
— Functionmat2vec(A::AbstractMatrix)
Converts a matrix to a vector.
sourceQuPhys.vec2mat
— Functionvec2mat(A::AbstractVector)
Converts a vector to a matrix.
sourceQuantum states, operators and super-operators
QuPhys.spre
— Functionspre(O::QuantumObject)
Returns the super-operator form of O
acting on the left of the density matrix operator, $\mathcal{O} \left(\hat{O}\right) \left[ \hat{\rho} \right] = \hat{O} \hat{\rho}$.
Since the density matrix is vectorized, this super-operator is always a matrix, obtained from $\mathcal{O} \left(\hat{O}\right) \boldsymbol{\cdot} = \hat{\mathbb{1}} \otimes \hat{O}$.
sourceQuPhys.spost
— Functionspost(O::QuantumObject)
Returns the super-operator form of O
acting on the right of the density matrix operator, $\mathcal{O} \left(\hat{O}\right) \left[ \hat{\rho} \right] = \hat{\rho} \hat{O}$.
Since the density matrix is vectorized, this super-operator is always a matrix, obtained from $\mathcal{O} \left(\hat{O}\right) \boldsymbol{\cdot} = \hat{O}^T \otimes \hat{\mathbb{1}}$.
sourceQuPhys.sprepost
— Functionsprepost(A::QuantumObject, B::QuantumObject)
Returns the super-operator form of A
and B
acting on the left and the right of the density matrix operator respectively, $\mathcal{O} \left( \hat{A}, \hat{B} \right) \left[ \hat{\rho} \right] = \hat{A} \hat{\rho} \hat{B}$.
Since the density matrix is vectorized, this super-operator is always a matrix, obtained from $\mathcal{O} \left(\hat{A}, \hat{B}\right) \boldsymbol{\cdot} = \text{spre}(A) * \text{spost}(B)$.
sourceQuPhys.lindblad_dissipator
— Functionlindblad_dissipator(O::QuantumObject)
Returns the Lindblad super-operator defined as $\mathcal{D} \left( \hat{O} \right) \left[ \hat{\rho} \right] = \frac{1}{2} \left( 2 \hat{O} \hat{\rho} \hat{O}^\dagger - \hat{O}^\dagger \hat{O} \hat{\rho} - \hat{\rho} \hat{O}^\dagger \hat{O} \right)$ considering the density matrix $\hat{\rho}$ in the vectorized form.
sourceQuPhys.destroy
— Functiondestroy(N::Int)
Bosonic annihilation operator with Hilbert space cutoff N
. This operator acts on a fock state as $\hat{a} \ket{n} = \sqrt{n} \ket{n-1}$.
Examples
julia> a = destroy(20)
+true
sourceQuPhys.wigner
— Functionwigner(state::QuantumObject, xvec::AbstractVector, yvec::AbstractVector; g::Real=√2,
+ solver::WignerSolver=WignerLaguerre())
Generates the Wigner quasipropability distribution of state
at points xvec + 1im * yvec
. The g
parameter is a scaling factor related to the value of $\hbar$ in the commutation relation $[x, y] = i \hbar$ via $\hbar=2/g^2$ giving the default value $\hbar=1$.
The solver
parameter can be either WignerLaguerre()
or WignerClenshaw()
. The former uses the Laguerre polynomial expansion of the Wigner function, while the latter uses the Clenshaw algorithm. The Laguerre expansion is faster for sparse matrices, while the Clenshaw algorithm is faster for dense matrices. The WignerLaguerre
solver has an optional parallel
parameter which defaults to true
and uses multithreading to speed up the calculation.
sourceQuPhys.get_coherence
— Functionget_coherence(ψ::QuantumObject)
Get the coherence value $\alpha$ by measuring the expectation value of the destruction operator $\hat{a}$ on the state $\ket{\psi}$.
It returns both $\alpha$ and the state $\ket{\delta_\psi} = \exp ( \bar{\alpha} \hat{a} - \alpha \hat{a}^\dagger )$. The latter corresponds to the quantum fulctuations around the coherent state $\ket{\alpha}$.
sourceQuPhys.n_th
— Functionn_th(ω::Number, T::Real)
Gives the mean number of excitations in a mode with frequency ω at temperature T: $n_{\rm th} (\omega, T) = \frac{1}{e^{\omega/T} - 1}$
sourceQuPhys.get_data
— Functionget_data(A::QuantumObject)
Returns the data of a QuantumObject.
sourceQuPhys.mat2vec
— Functionmat2vec(A::AbstractMatrix)
Converts a matrix to a vector.
sourceQuPhys.vec2mat
— Functionvec2mat(A::AbstractVector)
Converts a vector to a matrix.
sourceQuantum states, operators and super-operators
QuPhys.spre
— Functionspre(O::QuantumObject)
Returns the super-operator form of O
acting on the left of the density matrix operator, $\mathcal{O} \left(\hat{O}\right) \left[ \hat{\rho} \right] = \hat{O} \hat{\rho}$.
Since the density matrix is vectorized, this super-operator is always a matrix, obtained from $\mathcal{O} \left(\hat{O}\right) \boldsymbol{\cdot} = \hat{\mathbb{1}} \otimes \hat{O}$.
sourceQuPhys.spost
— Functionspost(O::QuantumObject)
Returns the super-operator form of O
acting on the right of the density matrix operator, $\mathcal{O} \left(\hat{O}\right) \left[ \hat{\rho} \right] = \hat{\rho} \hat{O}$.
Since the density matrix is vectorized, this super-operator is always a matrix, obtained from $\mathcal{O} \left(\hat{O}\right) \boldsymbol{\cdot} = \hat{O}^T \otimes \hat{\mathbb{1}}$.
sourceQuPhys.sprepost
— Functionsprepost(A::QuantumObject, B::QuantumObject)
Returns the super-operator form of A
and B
acting on the left and the right of the density matrix operator respectively, $\mathcal{O} \left( \hat{A}, \hat{B} \right) \left[ \hat{\rho} \right] = \hat{A} \hat{\rho} \hat{B}$.
Since the density matrix is vectorized, this super-operator is always a matrix, obtained from $\mathcal{O} \left(\hat{A}, \hat{B}\right) \boldsymbol{\cdot} = \text{spre}(A) * \text{spost}(B)$.
sourceQuPhys.lindblad_dissipator
— Functionlindblad_dissipator(O::QuantumObject)
Returns the Lindblad super-operator defined as $\mathcal{D} \left( \hat{O} \right) \left[ \hat{\rho} \right] = \frac{1}{2} \left( 2 \hat{O} \hat{\rho} \hat{O}^\dagger - \hat{O}^\dagger \hat{O} \hat{\rho} - \hat{\rho} \hat{O}^\dagger \hat{O} \right)$ considering the density matrix $\hat{\rho}$ in the vectorized form.
sourceQuPhys.destroy
— Functiondestroy(N::Int)
Bosonic annihilation operator with Hilbert space cutoff N
. This operator acts on a fock state as $\hat{a} \ket{n} = \sqrt{n} \ket{n-1}$.
Examples
julia> a = destroy(20)
Quantum Object: type=Operator dims=[20] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 19 stored entries:
⠈⠢⡀⠀⠀⠀⠀⠀⠀⠀
@@ -130,7 +130,7 @@
⠀⠀⠀⠀⠀⠀⠀⠀⠈⠢
julia> fock(20, 3)' * a * fock(20, 4)
-2.0 + 0.0im
sourceQuPhys.create
— Functioncreate(N::Int)
Bosonic creation operator with Hilbert space cutoff N
. This operator acts on a fock state as $\hat{a}^\dagger \ket{n} = \sqrt{n+1} \ket{n+1}$.
Examples
julia> a_d = create(20)
+2.0 + 0.0im
sourceQuPhys.create
— Functioncreate(N::Int)
Bosonic creation operator with Hilbert space cutoff N
. This operator acts on a fock state as $\hat{a}^\dagger \ket{n} = \sqrt{n+1} \ket{n+1}$.
Examples
julia> a_d = create(20)
Quantum Object: type=Operator dims=[20] size=(20, 20) ishermitian=false
20×20 SparseMatrixCSC{ComplexF64, Int64} with 19 stored entries:
⠢⡀⠀⠀⠀⠀⠀⠀⠀⠀
@@ -140,7 +140,7 @@
⠀⠀⠀⠀⠀⠀⠀⠈⠢⡀
julia> fock(20, 4)' * a_d * fock(20, 3)
-2.0 + 0.0im
sourceQuPhys.sigmap
— Functionsigmap()
Pauli ladder operator $\hat{\sigma}_+ = \hat{\sigma}_x + i \hat{\sigma}_y$.
sourceQuPhys.sigmam
— Functionsigmam()
Pauli ladder operator $\hat{\sigma}_- = \hat{\sigma}_x - i \hat{\sigma}_y$.
sourceQuPhys.sigmax
— Functionsigmax()
Pauli operator $\hat{\sigma}_x = \hat{\sigma}_- + \hat{\sigma}_+$.
sourceQuPhys.sigmay
— Functionsigmay()
Pauli operator $\hat{\sigma}_y = i \left( \hat{\sigma}_- - \hat{\sigma}_+ \right)$.
sourceQuPhys.sigmaz
— Functionsigmaz()
Pauli operator $\hat{\sigma}_z = \comm{\hat{\sigma}_+}{\hat{\sigma}_-}$.
sourceQuPhys.eye
— Functioneye(N::Int; type=OperatorQuantumObject, dims=[N])
Identity operator $\hat{\mathbb{1}}$ with Hilbert dimension N
.
sourceQuPhys.fock
— Functionfock(N::Int, pos::Int; dims::Vector{Int}=[N], sparse::Bool=false)
Generates a fock state $\ket{\psi}$ of dimension N
. It is also possible to specify the list of dimensions dims
if different subsystems are present.
sourceQuPhys.basis
— Functionbasis(N::Int, pos::Int; dims::Vector{Int}=[N])
Generates a fock state like fock
.
sourceQuPhys.coherent
— Functioncoherent(N::Real, α::T)
Generates a coherent state $\ket{\alpha}$, which is defined as an eigenvector of the bosonic annihilation operator $\hat{a} \ket{\alpha} = \alpha \ket{\alpha}$.
sourceQuPhys.rand_dm
— Functionrand_dm(N::Integer; kwargs...)
Generates a random density matrix $\hat{\rho}$, with the property to be positive definite, and that $\Tr \left[ \hat{\rho} \right] = 1$.
sourceQuPhys.projection
— Functionprojection(N::Int, i::Int, j::Int)
Generates the projection operator $\hat{O} = \dyad{i}{j}$ with Hilbert space dimension N
.
sourceQuPhys.sinm
— Functionsinm(O::QuantumObject)
Generates the sine of the operator O
, defined as
$\sin \left( \hat{O} \right) = \frac{e^{i \hat{O}} - e^{-i \hat{O}}}{2 i}$
sourceQuPhys.cosm
— Functioncosm(O::QuantumObject)
Generates the cosine of the operator O
, defined as
$\cos \left( \hat{O} \right) = \frac{e^{i \hat{O}} + e^{-i \hat{O}}}{2}$
sourceTime evolution
QuPhys.sesolveProblem
— FunctionsesolveProblem(H::QuantumObject,
+2.0 + 0.0im
sourceQuPhys.sigmap
— Functionsigmap()
Pauli ladder operator $\hat{\sigma}_+ = \hat{\sigma}_x + i \hat{\sigma}_y$.
sourceQuPhys.sigmam
— Functionsigmam()
Pauli ladder operator $\hat{\sigma}_- = \hat{\sigma}_x - i \hat{\sigma}_y$.
sourceQuPhys.sigmax
— Functionsigmax()
Pauli operator $\hat{\sigma}_x = \hat{\sigma}_- + \hat{\sigma}_+$.
sourceQuPhys.sigmay
— Functionsigmay()
Pauli operator $\hat{\sigma}_y = i \left( \hat{\sigma}_- - \hat{\sigma}_+ \right)$.
sourceQuPhys.sigmaz
— Functionsigmaz()
Pauli operator $\hat{\sigma}_z = \comm{\hat{\sigma}_+}{\hat{\sigma}_-}$.
sourceQuPhys.eye
— Functioneye(N::Int; type=OperatorQuantumObject, dims=[N])
Identity operator $\hat{\mathbb{1}}$ with Hilbert dimension N
.
sourceQuPhys.fock
— Functionfock(N::Int, pos::Int; dims::Vector{Int}=[N], sparse::Bool=false)
Generates a fock state $\ket{\psi}$ of dimension N
. It is also possible to specify the list of dimensions dims
if different subsystems are present.
sourceQuPhys.basis
— Functionbasis(N::Int, pos::Int; dims::Vector{Int}=[N])
Generates a fock state like fock
.
sourceQuPhys.coherent
— Functioncoherent(N::Real, α::T)
Generates a coherent state $\ket{\alpha}$, which is defined as an eigenvector of the bosonic annihilation operator $\hat{a} \ket{\alpha} = \alpha \ket{\alpha}$.
sourceQuPhys.rand_dm
— Functionrand_dm(N::Integer; kwargs...)
Generates a random density matrix $\hat{\rho}$, with the property to be positive definite, and that $\Tr \left[ \hat{\rho} \right] = 1$.
sourceQuPhys.projection
— Functionprojection(N::Int, i::Int, j::Int)
Generates the projection operator $\hat{O} = \dyad{i}{j}$ with Hilbert space dimension N
.
sourceQuPhys.sinm
— Functionsinm(O::QuantumObject)
Generates the sine of the operator O
, defined as
$\sin \left( \hat{O} \right) = \frac{e^{i \hat{O}} - e^{-i \hat{O}}}{2 i}$
sourceQuPhys.cosm
— Functioncosm(O::QuantumObject)
Generates the cosine of the operator O
, defined as
$\cos \left( \hat{O} \right) = \frac{e^{i \hat{O}} + e^{-i \hat{O}}}{2}$
sourceTime evolution
QuPhys.sesolveProblem
— FunctionsesolveProblem(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector;
alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm=Tsit5()
@@ -148,7 +148,7 @@
H_t::Union{Nothing,Function}=nothing,
params::Dict{Symbol, Any}=Dict{Symbol, Any}(),
progress::Bool=true,
- kwargs...)
Generates the ODEProblem for the Schrödinger time evolution of a quantum system.
Arguments
H::QuantumObject
: The Hamiltonian of the system.ψ0::QuantumObject
: The initial state of the system.t_l::AbstractVector
: The time list of the evolution.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: The algorithm used for the time evolution.e_ops::AbstractVector
: The list of operators to be evaluated during the evolution.H_t::Union{Nothing,Function}
: The time-dependent Hamiltonian of the system. If nothing
, the Hamiltonian is time-independent.params::Dict{Symbol, Any}
: The parameters of the system.progress::Bool
: Whether to show the progress bar.kwargs...
: The keyword arguments passed to the ODEProblem
constructor.
Returns
prob
: The ODEProblem
for the Schrödinger time evolution of the system.
sourceQuPhys.mesolveProblem
— FunctionmesolveProblem(H::QuantumObject,
+ kwargs...)
Generates the ODEProblem for the Schrödinger time evolution of a quantum system.
Arguments
H::QuantumObject
: The Hamiltonian of the system.ψ0::QuantumObject
: The initial state of the system.t_l::AbstractVector
: The time list of the evolution.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: The algorithm used for the time evolution.e_ops::AbstractVector
: The list of operators to be evaluated during the evolution.H_t::Union{Nothing,Function}
: The time-dependent Hamiltonian of the system. If nothing
, the Hamiltonian is time-independent.params::Dict{Symbol, Any}
: The parameters of the system.progress::Bool
: Whether to show the progress bar.kwargs...
: The keyword arguments passed to the ODEProblem
constructor.
Returns
prob
: The ODEProblem
for the Schrödinger time evolution of the system.
sourceQuPhys.mesolveProblem
— FunctionmesolveProblem(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector, c_ops::AbstractVector=[];
alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm=Tsit5(),
@@ -156,7 +156,7 @@
H_t::Union{Nothing,Function}=nothing,
params::Dict{Symbol, Any}=Dict{Symbol, Any}(),
progress::Bool=true,
- kwargs...)
Generates the ODEProblem for the master equation time evolution of an open quantum system.
Arguments
H::QuantumObject
: The Hamiltonian or the Liouvillian of the system.ψ0::QuantumObject
: The initial state of the system.t_l::AbstractVector
: The time list of the evolution.c_ops::AbstractVector=[]
: The list of the collapse operators.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm=Tsit5()
: The algorithm used for the time evolution.e_ops::AbstractVector=[]
: The list of the operators for which the expectation values are calculated.H_t::Union{Nothing,Function}=nothing
: The time-dependent Hamiltonian or Liouvillian.params::Dict{Symbol, Any}=Dict{Symbol, Any}()
: The parameters of the time evolution.progress::Bool=true
: Whether to show the progress bar.kwargs...
: The keyword arguments for the ODEProblem.
Returns
prob::ODEProblem
: The ODEProblem for the master equation time evolution.
sourceQuPhys.mcsolveProblem
— FunctionmcsolveProblem(H::QuantumObject,
+ kwargs...)
Generates the ODEProblem for the master equation time evolution of an open quantum system.
Arguments
H::QuantumObject
: The Hamiltonian or the Liouvillian of the system.ψ0::QuantumObject
: The initial state of the system.t_l::AbstractVector
: The time list of the evolution.c_ops::AbstractVector=[]
: The list of the collapse operators.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm=Tsit5()
: The algorithm used for the time evolution.e_ops::AbstractVector=[]
: The list of the operators for which the expectation values are calculated.H_t::Union{Nothing,Function}=nothing
: The time-dependent Hamiltonian or Liouvillian.params::Dict{Symbol, Any}=Dict{Symbol, Any}()
: The parameters of the time evolution.progress::Bool=true
: Whether to show the progress bar.kwargs...
: The keyword arguments for the ODEProblem.
Returns
prob::ODEProblem
: The ODEProblem for the master equation time evolution.
sourceQuPhys.mcsolveProblem
— FunctionmcsolveProblem(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector, c_ops::AbstractVector;
alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm=Tsit5(),
@@ -164,7 +164,7 @@
H_t::Union{Nothing,Function}=nothing,
params::Dict{Symbol, Any}=Dict{Symbol, Any}(),
jump_interp_pts::Integer=-1,
- kwargs...)
Generates the ODEProblem for a single trajectory of the Monte Carlo wave function time evolution of an open quantum system.
Arguments
H::QuantumObject
: Hamiltonian of the system.ψ0::QuantumObject
: Initial state of the system.t_l::AbstractVector
: List of times at which to save the state of the system.c_ops::AbstractVector
: List of collapse operators.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::AbstractVector
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.jump_interp_pts::Integer
: Number of points to use for interpolation of the jump times.kwargs...
: Additional keyword arguments to pass to the solver.
Returns
prob::ODEProblem
: The ODEProblem for the Monte Carlo wave function time evolution.
Notes
When jump_interp_pts
is set to -1, a DiscreteCallback
is used to detect the jump times. When jump_interp_pts
is set to a positive integer, a ContinuousCallback
is used to detect the jump times.
sourceQuPhys.mcsolveEnsembleProblem
— FunctionmcsolveEnsembleProblem(H::QuantumObject,
+ kwargs...)
Generates the ODEProblem for a single trajectory of the Monte Carlo wave function time evolution of an open quantum system.
Arguments
H::QuantumObject
: Hamiltonian of the system.ψ0::QuantumObject
: Initial state of the system.t_l::AbstractVector
: List of times at which to save the state of the system.c_ops::AbstractVector
: List of collapse operators.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::AbstractVector
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.jump_interp_pts::Integer
: Number of points to use for interpolation of the jump times.kwargs...
: Additional keyword arguments to pass to the solver.
Returns
prob::ODEProblem
: The ODEProblem for the Monte Carlo wave function time evolution.
Notes
When jump_interp_pts
is set to -1, a DiscreteCallback
is used to detect the jump times. When jump_interp_pts
is set to a positive integer, a ContinuousCallback
is used to detect the jump times.
sourceQuPhys.mcsolveEnsembleProblem
— FunctionmcsolveEnsembleProblem(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector, c_ops::AbstractVector;
alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm=Tsit5(),
@@ -176,7 +176,7 @@
jump_interp_pts::Integer=-1,
prob_func::Function=_mcsolve_prob_func,
output_func::Function=_mcsolve_output_func,
- kwargs...)
Generates the ODEProblem for an ensemble of trajectories of the Monte Carlo wave function time evolution of an open quantum system.
Arguments
H::QuantumObject
: Hamiltonian of the system.ψ0::QuantumObject
: Initial state of the system.t_l::AbstractVector
: List of times at which to save the state of the system.c_ops::AbstractVector
: List of collapse operators.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::AbstractVector
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.progress::Bool
: Whether to show a progress bar.n_traj::Integer
: Number of trajectories to use.jump_interp_pts::Integer
: Number of points to use for interpolation of the jump times.prob_func::Function
: Function to use for generating the ODEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.kwargs...
: Additional keyword arguments to pass to the solver.
Returns
prob::EnsembleProblem with ODEProblem
: The Ensemble ODEProblem for the Monte Carlo
wave function time evolution.
Notes
When jump_interp_pts
is set to -1, a DiscreteCallback
is used to detect the jump times. When jump_interp_pts
is set to a positive integer, a ContinuousCallback
is used to detect the jump times.
sourceQuPhys.sesolve
— Functionsesolve(H::QuantumObject,
+ kwargs...)
Generates the ODEProblem for an ensemble of trajectories of the Monte Carlo wave function time evolution of an open quantum system.
Arguments
H::QuantumObject
: Hamiltonian of the system.ψ0::QuantumObject
: Initial state of the system.t_l::AbstractVector
: List of times at which to save the state of the system.c_ops::AbstractVector
: List of collapse operators.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::AbstractVector
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.progress::Bool
: Whether to show a progress bar.n_traj::Integer
: Number of trajectories to use.jump_interp_pts::Integer
: Number of points to use for interpolation of the jump times.prob_func::Function
: Function to use for generating the ODEProblem.output_func::Function
: Function to use for generating the output of a single trajectory.kwargs...
: Additional keyword arguments to pass to the solver.
Returns
prob::EnsembleProblem with ODEProblem
: The Ensemble ODEProblem for the Monte Carlo
wave function time evolution.
Notes
When jump_interp_pts
is set to -1, a DiscreteCallback
is used to detect the jump times. When jump_interp_pts
is set to a positive integer, a ContinuousCallback
is used to detect the jump times.
sourceQuPhys.sesolve
— Functionsesolve(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector;
alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm=Tsit5(),
@@ -184,7 +184,7 @@
H_t::Union{Nothing,Function}=nothing,
params::Dict{Symbol, Any}=Dict{Symbol, Any}(),
progress::Bool=true,
- kwargs...)
Time evolution of a closed quantum system using the Schrödinger equation.
Arguments
H::QuantumObject
: Hamiltonian of the system.
ψ0::QuantumObject
: Initial state of the system.
t_l::AbstractVector
: List of times at which to save the state of the system.
alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.
e_ops::AbstractVector
: List of operators for which to calculate expectation values.
H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.
params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.
progress::Bool
: Whether to show a progress bar.
kwargs...
: Additional keyword arguments to pass to the solver.
Returns
sol::TimeEvolutionSol
: The solution of the time evolution.
sourceQuPhys.mesolve
— Functionmesolve(H::QuantumObject,
+ kwargs...)
Time evolution of a closed quantum system using the Schrödinger equation.
Arguments
H::QuantumObject
: Hamiltonian of the system.
ψ0::QuantumObject
: Initial state of the system.
t_l::AbstractVector
: List of times at which to save the state of the system.
alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.
e_ops::AbstractVector
: List of operators for which to calculate expectation values.
H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.
params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.
progress::Bool
: Whether to show a progress bar.
kwargs...
: Additional keyword arguments to pass to the solver.
Returns
sol::TimeEvolutionSol
: The solution of the time evolution.
sourceQuPhys.mesolve
— Functionmesolve(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector, c_ops::AbstractVector=[];
alg::OrdinaryDiffEqAlgorithm=Tsit5(),
@@ -192,7 +192,7 @@
H_t::Union{Nothing,Function}=nothing,
params::Dict{Symbol, Any}=Dict{Symbol, Any}(),
progress::Bool=true,
- kwargs...)
Time evolution of an open quantum system using master equation.
Arguments
H::QuantumObject
: Hamiltonian of Liouvillian of the system.ψ0::QuantumObject
: Initial state of the system.t_l::AbstractVector
: List of times at which to save the state of the system.c_ops::AbstractVector
: List of collapse operators.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::AbstractVector
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.progress::Bool
: Whether to show a progress bar.kwargs...
: Additional keyword arguments to pass to the solver.
Returns
sol::TimeEvolutionSol
: The solution of the time evolution.
sourceQuPhys.mcsolve
— Functionmcsolve(H::QuantumObject,
+ kwargs...)
Time evolution of an open quantum system using master equation.
Arguments
H::QuantumObject
: Hamiltonian of Liouvillian of the system.ψ0::QuantumObject
: Initial state of the system.t_l::AbstractVector
: List of times at which to save the state of the system.c_ops::AbstractVector
: List of collapse operators.alg::OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::AbstractVector
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.progress::Bool
: Whether to show a progress bar.kwargs...
: Additional keyword arguments to pass to the solver.
Returns
sol::TimeEvolutionSol
: The solution of the time evolution.
sourceQuPhys.mcsolve
— Functionmcsolve(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector, c_ops::AbstractVector;
alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm=Tsit5(),
@@ -203,7 +203,7 @@
n_traj::Integer=1,
ensemble_method=EnsembleThreads(),
jump_interp_pts::Integer=-1,
- kwargs...)
Time evolution of an open quantum system using quantum trajectories.
Arguments
H::QuantumObject
: Hamiltonian of the system.ψ0::QuantumObject
: Initial state of the system.t_l::AbstractVector
: List of times at which to save the state of the system.c_ops::AbstractVector
: List of collapse operators.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::AbstractVector
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.progress::Bool
: Whether to show a progress bar.n_traj::Integer
: Number of trajectories to use.ensemble_method
: Ensemble method to use.jump_interp_pts::Integer
: Number of points to use for interpolation of jump times.kwargs...
: Additional keyword arguments to pass to the solver.
Returns
sol::TimeEvolutionMCSol
: The solution of the time evolution.
Notes
ensemble_method
can be one of EnsembleThreads()
, EnsembleSerial()
, EnsembleDistributed()
. When jump_interp_pts
is set to -1, a DiscreteCallback
is used to detect the jump times. When jump_interp_pts
is set to a positive integer, a ContinuousCallback
is used to detect the jump times.
sourceQuPhys.dfd_mesolve
— Functionfunction dfd_mesolve(H::Function, ψ0::QuantumObject,
+ kwargs...)
Time evolution of an open quantum system using quantum trajectories.
Arguments
H::QuantumObject
: Hamiltonian of the system.ψ0::QuantumObject
: Initial state of the system.t_l::AbstractVector
: List of times at which to save the state of the system.c_ops::AbstractVector
: List of collapse operators.alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm
: Algorithm to use for the time evolution.e_ops::AbstractVector
: List of operators for which to calculate expectation values.H_t::Union{Nothing,Function}
: Time-dependent part of the Hamiltonian.params::Dict{Symbol, Any}
: Dictionary of parameters to pass to the solver.progress::Bool
: Whether to show a progress bar.n_traj::Integer
: Number of trajectories to use.ensemble_method
: Ensemble method to use.jump_interp_pts::Integer
: Number of points to use for interpolation of jump times.kwargs...
: Additional keyword arguments to pass to the solver.
Returns
sol::TimeEvolutionMCSol
: The solution of the time evolution.
Notes
ensemble_method
can be one of EnsembleThreads()
, EnsembleSerial()
, EnsembleDistributed()
. When jump_interp_pts
is set to -1, a DiscreteCallback
is used to detect the jump times. When jump_interp_pts
is set to a positive integer, a ContinuousCallback
is used to detect the jump times.
sourceQuPhys.dfd_mesolve
— Functionfunction dfd_mesolve(H::Function, ψ0::QuantumObject,
t_l::AbstractVector, c_ops::Function, maxdims::AbstractVector;
alg::OrdinaryDiffEq.OrdinaryDiffEqAlgorithm=Tsit5(),
e_ops::Union{Nothing, Function}=nothing,
@@ -211,7 +211,7 @@
params::Dict{Symbol, Any}=Dict{Symbol, Any}(),
progress::Bool=true,
tol_list::Vector{<:Number}=Float64[],
- kwargs...)
Time evolution of an open quantum system using master equation, dynamically changing the dimension of the Hilbert subspaces.
sourceQuPhys.dsf_mesolve
— Functionfunction dsf_mesolve(H::Function,
+ kwargs...)
Time evolution of an open quantum system using master equation, dynamically changing the dimension of the Hilbert subspaces.
sourceQuPhys.dsf_mesolve
— Functionfunction dsf_mesolve(H::Function,
ψ0::QuantumObject,
t_l::AbstractVector, c_ops::Function,
op_list::AbstractVector,
@@ -222,7 +222,7 @@
params::Dict{Symbol, Any}=Dict{Symbol, Any}(),
progress::Bool=true,
δα_list::Vector{<:Number}=Float64[],
- kwargs...)
Time evolution of an open quantum system using master equation and the Dynamical Shifted Fock algorithm.
sourceQuPhys.dsf_mcsolve
— Functionfunction dsf_mcsolve(H::Function,
+ kwargs...)
Time evolution of an open quantum system using master equation and the Dynamical Shifted Fock algorithm.
sourceQuPhys.dsf_mcsolve
— Functionfunction dsf_mcsolve(H::Function,
ψ0::QuantumObject,
t_l::AbstractVector, c_ops::Function,
op_list::AbstractVector,
@@ -237,12 +237,12 @@
ensemble_method=EnsembleThreads(),
jump_interp_pts::Integer=10,
krylov_dim::Integer=cld(prod(ψ0.dims), 4),
- kwargs...)
Time evolution of a quantum system using the Monte Carlo wave function method and the Dynamical Shifted Fock algorithm.
sourceQuPhys.liouvillian
— Functionliouvillian(H::QuantumObject, c_ops::AbstractVector)
Construct the Liouvillian superoperator for a system Hamiltonian and a set of collapse operators: $\mathcal{L} \cdot = -i[\hat{H}, \cdot] + \sum_i \mathcal{D}[\hat{O}_i] \cdot$, where $\mathcal{D}[\hat{O}_i] \cdot = \hat{O}_i \cdot \hat{O}_i^\dagger - \frac{1}{2} \hat{O}_i^\dagger \hat{O}_i \cdot - \frac{1}{2} \cdot \hat{O}_i^\dagger \hat{O}_i$.
sourceQuPhys.liouvillian_generalized
— Functionliouvillian_generalized(H::QuantumObject, fields::Vector,
-T_list::Vector; N_trunc::Int=size(H,1), tol::Float64=0.0, σ_filter::Union{Nothing, Real}=nothing)
Constructs the generalized Liouvillian for a system coupled to a bath of harmonic oscillators.
See, e.g., Settineri, Alessio, et al. "Dissipation and thermal noise in hybrid quantum systems in the ultrastrong-coupling regime." Physical Review A 98.5 (2018): 053834.
sourceQuPhys.steadystate_floquet
— Functionsteadystate_floquet(H_0::QuantumObject,
+ kwargs...)
Time evolution of a quantum system using the Monte Carlo wave function method and the Dynamical Shifted Fock algorithm.
sourceQuPhys.liouvillian
— Functionliouvillian(H::QuantumObject, c_ops::AbstractVector)
Construct the Liouvillian superoperator for a system Hamiltonian and a set of collapse operators: $\mathcal{L} \cdot = -i[\hat{H}, \cdot] + \sum_i \mathcal{D}[\hat{O}_i] \cdot$, where $\mathcal{D}[\hat{O}_i] \cdot = \hat{O}_i \cdot \hat{O}_i^\dagger - \frac{1}{2} \hat{O}_i^\dagger \hat{O}_i \cdot - \frac{1}{2} \cdot \hat{O}_i^\dagger \hat{O}_i$.
sourceQuPhys.liouvillian_generalized
— Functionliouvillian_generalized(H::QuantumObject, fields::Vector,
+T_list::Vector; N_trunc::Int=size(H,1), tol::Float64=0.0, σ_filter::Union{Nothing, Real}=nothing)
Constructs the generalized Liouvillian for a system coupled to a bath of harmonic oscillators.
See, e.g., Settineri, Alessio, et al. "Dissipation and thermal noise in hybrid quantum systems in the ultrastrong-coupling regime." Physical Review A 98.5 (2018): 053834.
sourceQuPhys.steadystate_floquet
— Functionsteadystate_floquet(H_0::QuantumObject,
c_ops::Vector, H_p::QuantumObject,
H_m::QuantumObject,
ω::Real; n_max::Int=4, lf_solver::LSolver=LiouvillianDirectSolver(),
- ss_solver::Type{SSSolver}=SteadyStateDirectSolver)
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 $\\omega$, we divide it into two parts, H_p
and H_m
, where H_p
oscillates as $e^{i \\omega t}$ and H_m
oscillates as $e^{-i \\omega t}$. n_max
is the number of iterations used to obtain the effective Liouvillian, lf_solver
is the solver used to solve the effective Liouvillian, and ss_solver
is the solver used to solve the steady state.
sourceCorrelations and Spectrum
QuPhys.correlation_3op_2t
— Functioncorrelation_3op_2t(H::QuantumObject,
+ ss_solver::Type{SSSolver}=SteadyStateDirectSolver)
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 $\\omega$, we divide it into two parts, H_p
and H_m
, where H_p
oscillates as $e^{i \\omega t}$ and H_m
oscillates as $e^{-i \\omega t}$. n_max
is the number of iterations used to obtain the effective Liouvillian, lf_solver
is the solver used to solve the effective Liouvillian, and ss_solver
is the solver used to solve the steady state.
sourceCorrelations and Spectrum
QuPhys.correlation_3op_2t
— Functioncorrelation_3op_2t(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector,
τ_l::AbstractVector,
@@ -250,7 +250,7 @@
B::QuantumObject,
C::QuantumObject,
c_ops::AbstractVector=[];
- kwargs...)
Returns the two-times correlation function of three operators $\hat{A}$, $\hat{B}$ and $\hat{C}$: $\expval{\hat{A}(t) \hat{B}(t + \tau) \hat{C}(t)}$ for a given initial state $\ket{\psi_0}$.
sourceQuPhys.correlation_2op_2t
— Functioncorrelation_2op_2t(H::QuantumObject,
+ kwargs...)
Returns the two-times correlation function of three operators $\hat{A}$, $\hat{B}$ and $\hat{C}$: $\expval{\hat{A}(t) \hat{B}(t + \tau) \hat{C}(t)}$ for a given initial state $\ket{\psi_0}$.
sourceQuPhys.correlation_2op_2t
— Functioncorrelation_2op_2t(H::QuantumObject,
ψ0::QuantumObject,
t_l::AbstractVector,
τ_l::AbstractVector,
@@ -258,20 +258,20 @@
B::QuantumObject,
c_ops::AbstractVector=[];
reverse::Bool=false,
- kwargs...)
Returns the two-times correlation function of two operators $\hat{A}$ and $\hat{B}$at different times
\expval{\hat{A}(t + \tau) \hat{B}(t)}
. When
reverse=true
, the correlation function is calculated as
\expval{\hat{A}(t) \hat{B}(t + \tau)}
`.
sourceQuPhys.correlation_2op_1t
— Functioncorrelation_2op_1t(H::QuantumObject,
+ kwargs...)
Returns the two-times correlation function of two operators $\hat{A}$ and $\hat{B}$at different times
\expval{\hat{A}(t + \tau) \hat{B}(t)}
. When
reverse=true
, the correlation function is calculated as
\expval{\hat{A}(t) \hat{B}(t + \tau)}
`.
sourceQuPhys.correlation_2op_1t
— Functioncorrelation_2op_1t(H::QuantumObject,
ψ0::QuantumObject,
τ_l::AbstractVector,
A::QuantumObject,
B::QuantumObject,
c_ops::AbstractVector=[];
reverse::Bool=false,
- kwargs...)
Returns the one-time correlation function of two operators $\hat{A}$ and $\hat{B}$ at different times $\expval{\hat{A}(\tau) \hat{B}(0)}$. When $reverse=true$, the correlation function is calculated as $\expval{\hat{A}(0) \hat{B}(\tau)}$.
sourceQuPhys.spectrum
— Functionspectrum(H::QuantumObject,
+ kwargs...)
Returns the one-time correlation function of two operators $\hat{A}$ and $\hat{B}$ at different times $\expval{\hat{A}(\tau) \hat{B}(0)}$. When $reverse=true$, the correlation function is calculated as $\expval{\hat{A}(0) \hat{B}(\tau)}$.
sourceQuPhys.spectrum
— Functionspectrum(H::QuantumObject,
ω_list::AbstractVector,
A::QuantumObject{<:AbstractArray{T2},OperatorQuantumObject},
B::QuantumObject{<:AbstractArray{T3},OperatorQuantumObject},
c_ops::AbstractVector=[];
solver::MySolver=ExponentialSeries(),
- kwargs...)
Returns the emission spectrum $S(\omega) = \int_{-\infty}^\infty \expval{\hat{A}(\tau) \hat{B}(0)} e^{-i \omega \tau} d \tau$.
sourceEigenvalues and eigenvectors
LinearAlgebra.eigen
— FunctionLinearAlgebra.eigen(A::QuantumObject; kwargs...)
Calculates the eigenvalues and eigenvectors of the QuantumObject
A
using the Julia LinearAlgebra package.
julia> a = destroy(5);
+ kwargs...)
Returns the emission spectrum $S(\omega) = \int_{-\infty}^\infty \expval{\hat{A}(\tau) \hat{B}(0)} e^{-i \omega \tau} d \tau$.
sourceEigenvalues and eigenvectors
LinearAlgebra.eigen
— FunctionLinearAlgebra.eigen(A::QuantumObject; kwargs...)
Calculates the eigenvalues and eigenvectors of the QuantumObject
A
using the Julia LinearAlgebra package.
julia> a = destroy(5);
julia> H = a + a'
Quantum Object: type=Operator dims=[5] size=(5, 5) ishermitian=true
@@ -302,10 +302,10 @@
julia> ψ_1 = QuantumObject(U[:,1], dims=H.dims);
julia> expect(H, ψ_1) ≈ E[1]
-true
sourceLinearAlgebra.eigvals
— FunctionLinearAlgebra.eigvals(A::QuantumObject; kwargs...)
Same as eigen(A::QuantumObject; kwargs...)
but for only the eigenvalues.
sourceQuPhys.eigsolve
— Functionfunction eigsolve(A::QuantumObject; v0::Union{Nothing,AbstractVector}=nothing,
+true
sourceLinearAlgebra.eigvals
— FunctionLinearAlgebra.eigvals(A::QuantumObject; kwargs...)
Same as eigen(A::QuantumObject; kwargs...)
but for only the eigenvalues.
sourceQuPhys.eigsolve
— Functionfunction eigsolve(A::QuantumObject; v0::Union{Nothing,AbstractVector}=nothing,
sigma::Union{Nothing, Real}=nothing, k::Integer = min(4, size(A, 1)),
krylovdim::Integer = min(10, size(A, 1)), tol::Real = 1e-8, maxiter::Integer = 200,
- solver::Union{Nothing, LinearSolve.SciMLLinearSolveAlgorithm} = nothing, showinfo::Bool=false, kwargs...)
Solve for the eigenvalues and eigenvectors of a matrix A
using the Arnoldi method. The keyword arguments are passed to the linear solver.
sourceQuPhys.eigsolve_al
— Functioneigsolve_al(H::QuantumObject,
+ solver::Union{Nothing, LinearSolve.SciMLLinearSolveAlgorithm} = nothing, showinfo::Bool=false, kwargs...)
Solve for the eigenvalues and eigenvectors of a matrix A
using the Arnoldi method. The keyword arguments are passed to the linear solver.
sourceQuPhys.eigsolve_al
— Functioneigsolve_al(H::QuantumObject,
T::Real, c_ops::AbstractVector=[];
alg::OrdinaryDiffEqAlgorithm=Tsit5(),
H_t::Union{Nothing,Function}=nothing,
@@ -317,4 +317,4 @@
maxiter::Integer=200,
eigstol::Real=1e-6,
showinfo::Bool=false,
- kwargs...)
Solve the eigenvalue problem for a Liouvillian superoperator L
using the Arnoldi-Lindblad method.
Arguments
H
: The Hamiltonian (or directly the Liouvillian) of the system.T
: The time at which to evaluate the time evolutionc_ops
: A vector of collapse operatorsalg
: The differential equation solver algorithmH_t
: A function H_t(t)
that returns the additional term at time t
params
: A dictionary of additional parametersprogress
: Whether to show a progress barρ0
: The initial density matrix. If not specified, a random density matrix is usedk
: The number of eigenvalues to computekrylovdim
: The dimension of the Krylov subspacemaxiter
: The maximum number of iterations for the eigsolvereigstol
: The tolerance for the eigsolvershowinfo
: Whether to show information about the eigsolverkwargs
: Additional keyword arguments passed to the differential equation solver
Returns
vals
: The eigenvaluesvecs
: The eigenvectors
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.
sourceSettings
This document was generated with Documenter.jl version 1.1.0 on Sunday 8 October 2023. Using Julia version 1.9.3.
+ kwargs...)
Solve the eigenvalue problem for a Liouvillian superoperator L
using the Arnoldi-Lindblad method.
Arguments
H
: The Hamiltonian (or directly the Liouvillian) of the system.T
: The time at which to evaluate the time evolutionc_ops
: A vector of collapse operatorsalg
: The differential equation solver algorithmH_t
: A functionH_t(t)
that returns the additional term at timet
params
: A dictionary of additional parametersprogress
: Whether to show a progress barρ0
: The initial density matrix. If not specified, a random density matrix is usedk
: The number of eigenvalues to computekrylovdim
: The dimension of the Krylov subspacemaxiter
: The maximum number of iterations for the eigsolvereigstol
: The tolerance for the eigsolvershowinfo
: Whether to show information about the eigsolverkwargs
: Additional keyword arguments passed to the differential equation solver
Returns
vals
: The eigenvaluesvecs
: The eigenvectors
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.