-
Notifications
You must be signed in to change notification settings - Fork 16
/
general_functions.jl
381 lines (306 loc) · 11.8 KB
/
general_functions.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
"""
row_major_reshape(Q::AbstractArray, shapes...)
Reshapes `Q` in the row-major order, as numpy.
"""
row_major_reshape(Q::AbstractArray{T}, shapes...) where {T} = PermutedDimsArray(reshape(Q, reverse(shapes)...), (length(shapes):-1:1))
"""
meshgrid(x::AbstractVector, y::AbstractVector)
Equivalent to [numpy meshgrid](https://numpy.org/doc/stable/reference/generated/numpy.meshgrid.html).
"""
function meshgrid(x::AbstractVector{T}, y::AbstractVector{T}) where {T}
X = reshape(repeat(x, inner=length(y)), length(y), length(x))
Y = repeat(y, outer=(1, length(x)))
X, Y
end
"""
sparse_to_dense(A::QuantumObject)
Converts a sparse QuantumObject to a dense QuantumObject.
"""
sparse_to_dense(A::QuantumObject{<:AbstractVecOrMat}) = QuantumObject(sparse_to_dense(A.data), A.type, A.dims)
sparse_to_dense(A::AbstractSparseArray) = Array(A)
for op in (:Transpose, :Adjoint)
@eval sparse_to_dense(A::$op{T, <:AbstractSparseMatrix}) where {T<:BlasFloat} = Array(A)
end
sparse_to_dense(A::AbstractArray) = A
"""
dense_to_sparse(A::QuantumObject)
Converts a dense QuantumObject to a sparse QuantumObject.
"""
dense_to_sparse(A::QuantumObject{<:AbstractVecOrMat}, tol::Real=1e-10) = QuantumObject(dense_to_sparse(A.data, tol), A.type, A.dims)
function dense_to_sparse(A::AbstractMatrix, tol::Real=1e-10)
idxs = findall(@. abs(A) > tol)
row_indices = getindex.(idxs, 1)
col_indices = getindex.(idxs, 2)
vals = getindex(A, idxs)
return sparse(row_indices, col_indices, vals, size(A)...)
end
function dense_to_sparse(A::AbstractVector, tol::Real=1e-10)
idxs = findall(@. abs(A) > tol)
vals = getindex(A, idxs)
return sparsevec(idxs, vals, length(A))
end
"""
tidyup(A::QuantumObject, tol::Real=1e-14)
Removes those elements of a QuantumObject `A` whose absolute value is less than `tol`.
"""
tidyup(A::QuantumObject{<:AbstractArray{T}}, tol::T2=1e-14) where {T, T2<:Real} = QuantumObject(tidyup(A.data, tol), A.type, A.dims)
tidyup(A::AbstractArray{T}, tol::T2=1e-14) where {T, T2<:Real} = @. T(abs(A) > tol) * A
tidyup(A::AbstractSparseMatrix{T}, tol::T2=1e-14) where {T, T2<:Real} = droptol!(copy(A), tol)
"""
tidyup!(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`](#tidyup).
"""
tidyup!(A::QuantumObject{<:AbstractArray{T}}, tol::T2=1e-14) where {T, T2<:Real} = (tidyup!(A.data, tol); A)
tidyup!(A::AbstractArray{T}, tol::T2=1e-14) where {T, T2<:Real} = @. A = T(abs(A) > tol) * A
tidyup!(A::AbstractSparseMatrix{T}, tol::T2=1e-14) where {T, T2<:Real} = droptol!(A, tol)
"""
get_data(A::QuantumObject)
Returns the data of a QuantumObject.
"""
get_data(A::QuantumObject) = A.data
"""
mat2vec(A::AbstractMatrix)
Converts a matrix to a vector.
"""
mat2vec(A::AbstractMatrix) = vec(A) # reshape(A, :)
function mat2vec(A::AbstractSparseMatrix)
i, j, v = findnz(A)
return sparsevec(i .+ (j .- 1) .* size(A, 1), v, prod(size(A)))
end
for op in (:Transpose, :Adjoint)
@eval mat2vec(A::$op{T, <:AbstractSparseMatrix}) where {T<:BlasFloat} = mat2vec(sparse(A))
@eval mat2vec(A::$op{T, <:AbstractMatrix}) where {T<:BlasFloat} = mat2vec(Matrix(A))
end
"""
vec2mat(A::AbstractVector)
Converts a vector to a matrix.
"""
function vec2mat(A::AbstractVector)
newsize = isqrt(length(A))
reshape(A, newsize, newsize)
end
@doc raw"""
gaussian(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.
"""
gaussian(x::Number, μ::Number, σ::Number) = exp(-(x - μ)^2 / (2 * σ^2))
@doc raw"""
ptrace(QO::QuantumObject, sel::Vector{Int})
[Partial trace](https://en.wikipedia.org/wiki/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
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 ``\ket{\psi} = \frac{1}{\sqrt{2}} \left( \ket{e,e} + \ket{g,g} \right)``:
```
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
```
"""
function ptrace(QO::QuantumObject{<:AbstractArray{T1},OpType}, sel::Vector{T2}) where
{T1,T2<:Int,OpType<:Union{KetQuantumObject,OperatorQuantumObject}}
length(QO.dims) == 1 && return QO
if isket(QO) || isbra(QO)
ρtr, dkeep = _ptrace_ket(QO.data, QO.dims, sel)
return QuantumObject(ρtr, dims=dkeep)
elseif isoper(QO)
ρtr, dkeep = _ptrace_oper(QO.data, QO.dims, sel)
return QuantumObject(ρtr, dims=dkeep)
end
end
@doc raw"""
entropy_vn(ρ::QuantumObject; base::Int=0, tol::Real=1e-15)
Calculates the [Von Neumann entropy](https://en.wikipedia.org/wiki/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
0.0 + 0.0im
julia> ρ = ket2dm(ψ)
Quantum Object: type=Operator dims=[2] size=(2, 2) ishermitian=true
2×2 Matrix{ComplexF64}:
1.0+0.0im 0.0+0.0im
0.0+0.0im 0.0+0.0im
julia> entropy_vn(ρ, base=2)
-0.0
```
Mixed state:
```
julia> ρ = 1 / 2 * ( ket2dm(fock(2,0)) + ket2dm(fock(2,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
julia> entropy_vn(ρ, base=2)
1.0
```
"""
function entropy_vn(ρ::QuantumObject{<:AbstractArray{T},OperatorQuantumObject}; base::Int=0, tol::Real=1e-15) where {T}
vals = eigvals(ρ)
indexes = abs.(vals) .> tol
1 ∉ indexes && return 0
nzvals = vals[indexes]
logvals = base != 0 ? log.(base, Complex.(nzvals)) : log.(Complex.(nzvals))
return -real(sum(nzvals .* logvals))
end
"""
entanglement(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`](@ref).
"""
function entanglement(QO::QuantumObject{<:AbstractArray{T},OpType}, sel::Vector{Int}) where
{T,OpType<:Union{BraQuantumObject,KetQuantumObject,OperatorQuantumObject}}
ψ = normalize(QO)
ρ_tr = ptrace(ψ, sel)
entropy = entropy_vn(ρ_tr)
return (entropy > 0) * entropy
end
@doc raw"""
expect(O::QuantumObject, ψ::QuantumObject)
Expectation value of the operator `O` with the state `ψ`. The latter
can be a [`KetQuantumObject`](@ref), [`BraQuantumObject`](@ref) or a [`OperatorQuantumObject`](@ref).
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
```
"""
function expect(O::QuantumObject{<:AbstractArray{T1},OperatorQuantumObject}, ψ::QuantumObject{<:AbstractArray{T2},KetQuantumObject}) where {T1,T2}
ψd = ψ.data
Od = O.data
return ishermitian(O) ? real(dot(ψd, Od, ψd)) : dot(ψd, Od, ψd)
end
function expect(O::QuantumObject{<:AbstractArray{T1},OperatorQuantumObject}, ψ::QuantumObject{<:AbstractArray{T2},BraQuantumObject}) where {T1,T2}
ψd = ψ.data'
Od = O.data
return ishermitian(O) ? real(dot(ψd, Od, ψd)) : dot(ψd, Od, ψd)
end
function expect(O::QuantumObject{<:AbstractArray{T1},OperatorQuantumObject}, ρ::QuantumObject{<:AbstractArray{T2},OperatorQuantumObject}) where {T1,T2}
return ishermitian(O) ? real(tr(O * ρ)) : tr(O * ρ)
end
@doc raw"""
get_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}``.
"""
function get_coherence(ψ::QuantumObject{<:AbstractArray{T}, StateOpType}) where {T,StateOpType<:Union{KetQuantumObject,OperatorQuantumObject}}
a = destroy(size(ψ,1))
α = expect(a, ψ)
D = exp(α*a' - conj(α)*a)
α, D' * ψ
end
@doc raw"""
n_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}``
"""
function n_th(ω::Real, T::Real)::Float64
(T == 0 || ω == 0) && return 0.0
abs(ω / T) > 50 && return 0.0
return 1 / (exp(ω / T) - 1)
end
function _ptrace_ket(QO::AbstractArray{T1}, dims::Vector{<:Integer}, sel::Vector{T2}) where
{T1,T2<:Integer}
rd = dims
nd = length(rd)
nd == 1 && return QO, rd
dkeep = rd[sel]
qtrace = filter!(e -> e ∉ sel, Vector(1:nd))
dtrace = @view(rd[qtrace])
vmat = reshape(QO, reverse(rd)...)
topermute = nd+1 .- vcat(sel, qtrace)
vmat = PermutedDimsArray(vmat, topermute)
vmat = reshape(vmat, prod(dkeep), prod(dtrace))
return vmat * vmat', dkeep
end
function _ptrace_oper(QO::AbstractArray{T1}, dims::Vector{<:Integer}, sel::Vector{T2}) where
{T1,T2<:Integer}
rd = dims
nd = length(rd)
nd == 1 && return QO, rd
dkeep = rd[sel]
qtrace = filter!(e -> e ∉ sel, Vector(1:nd))
dtrace = @view(rd[qtrace])
ρmat = reshape(QO, reverse!(repeat(rd, 2))...)
topermute = 2*nd+1 .- vcat(qtrace, qtrace .+ nd, sel, sel .+ nd)
reverse!(topermute)
ρmat = PermutedDimsArray(ρmat, topermute)
ρmat = row_major_reshape(ρmat, prod(dtrace), prod(dtrace), prod(dkeep), prod(dkeep))
res = dropdims(mapslices(tr, ρmat, dims=(1,2)), dims=(1,2))
return res, dkeep
end
function sparse_to_dense(::Type{M}) where M <: SparseMatrixCSC
T = M
par = T.parameters
npar = length(par)
(2 == npar) || error("Type $M is not supported.")
return Matrix{par[1]}
end
function mat2vec(::Type{M}) where M <: DenseMatrix
T = hasproperty(M, :body) ? M.body : M
par = T.parameters
npar = length(par)
(2 ≤ npar ≤ 3) || error("Type $M is not supported.")
if npar == 2
S = T.name.wrapper{par[1], 1}
else
S = T.name.wrapper{par[1], 1, par[3]}
end
return S
end
function mat2vec(::Type{M}) where M <: SparseMatrixCSC
T = M
par = T.parameters
npar = length(par)
(2 == npar) || error("Type $M is not supported.")
return SparseVector{par[1], par[2]}
end
function mat2vec(::Type{M}) where M <: Union{Adjoint{<:BlasFloat,<:SparseMatrixCSC}, Transpose{<:BlasFloat,<:SparseMatrixCSC}}
T = M.parameters[2]
par = T.parameters
npar = length(par)
(2 == npar) || error("Type $M is not supported.")
return SparseVector{par[1], par[2]}
end