Skip to content

Commit 7cac2f8

Browse files
authored
Make sure state generating functions return dense array (#591)
1 parent deb5d49 commit 7cac2f8

File tree

3 files changed

+15
-10
lines changed

3 files changed

+15
-10
lines changed

CHANGELOG.md

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -11,6 +11,7 @@ and this project adheres to [Semantic Versioning](https://semver.org/spec/v2.0.0
1111
- Add keyword argument `assume_hermitian` to `liouvillian`. This allows users to disable the assumption that the Hamiltonian is Hermitian. ([#581])
1212
- Use LinearSolve's internal methods for preconditioners in `SteadyStateLinearSolver`. ([#588])
1313
- Use `FillArrays.jl` for handling superoperators. This makes the code cleaner and potentially more efficient. ([#589])
14+
- Make sure state generating functions return dense array by default. ([#591])
1415

1516
## [v0.38.1]
1617
Release date: 2025-10-27
@@ -366,3 +367,4 @@ Release date: 2024-11-13
366367
[#581]: https://github.com/qutip/QuantumToolbox.jl/issues/581
367368
[#588]: https://github.com/qutip/QuantumToolbox.jl/issues/588
368369
[#589]: https://github.com/qutip/QuantumToolbox.jl/issues/589
370+
[#591]: https://github.com/qutip/QuantumToolbox.jl/issues/591

src/entropy.jl

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -36,9 +36,9 @@ Mixed state:
3636
julia> ρ = maximally_mixed_dm(2)
3737
3838
Quantum Object: type=Operator() dims=[2] size=(2, 2) ishermitian=true
39-
2×2 Diagonal{ComplexF64, Vector{ComplexF64}}:
40-
0.5-0.0im
41-
0.5-0.0im
39+
2×2 Matrix{ComplexF64}:
40+
0.5+0.0im 0.0+0.0im
41+
0.0+0.0im 0.5+0.0im
4242
4343
julia> entropy_vn(ρ, base=2)
4444
1.0

src/qobj/states.jl

Lines changed: 10 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -127,8 +127,7 @@ Density matrix for a thermal state (generating thermal state probabilities) with
127127
"""
128128
function thermal_dm(N::Int, n::Real; sparse::Union{Bool,Val} = Val(false))
129129
β = log(1.0 / n + 1.0)
130-
N_list = Array{Float64}(0:(N-1))
131-
data = exp.(-β .* N_list)
130+
data = [exp(-β * ComplexF64(j)) for j in 0:(N-1)]
132131
if getVal(sparse)
133132
return QuantumObject(spdiagm(0 => data ./ sum(data)), Operator(), N)
134133
else
@@ -149,10 +148,10 @@ The `dimensions` can be either the following types:
149148
If you want to keep type stability, it is recommended to use `maximally_mixed_dm(dimensions)` with `dimensions` as `Tuple` or `SVector` from [StaticArrays.jl](https://github.com/JuliaArrays/StaticArrays.jl) to keep type stability. See the [related Section](@ref doc:Type-Stability) about type stability for more details.
150149
"""
151150
maximally_mixed_dm(dimensions::Int) =
152-
QuantumObject(Eye(dimensions) / complex(dimensions), Operator(), SVector(dimensions))
151+
QuantumObject(diagm(0 => fill(ComplexF64(1 / dimensions), dimensions)), Operator(), SVector(dimensions))
153152
function maximally_mixed_dm(dimensions::Union{Dimensions,AbstractVector{Int},Tuple})
154153
N = prod(dimensions)
155-
return QuantumObject(Eye(N) / complex(N), Operator(), dimensions)
154+
return QuantumObject(diagm(0 => fill(ComplexF64(1 / N), N)), Operator(), dimensions)
156155
end
157156

158157
@doc raw"""
@@ -321,7 +320,9 @@ Returns the `n`-qubit [W-state](https://en.wikipedia.org/wiki/W_state):
321320
function w_state(::Val{n}) where {n}
322321
nzind = 2 .^ (0:(n-1)) .+ 1
323322
nzval = fill(ComplexF64(1 / sqrt(n)), n)
324-
return QuantumObject(SparseVector(2^n, nzind, nzval), Ket(), ntuple(x -> 2, Val(n)))
323+
data = zeros(ComplexF64, 2^n)
324+
@inbounds data[nzind] .= nzval
325+
return QuantumObject(data, Ket(), ntuple(x -> 2, Val(n)))
325326
end
326327
w_state(n::Int) = w_state(Val(n))
327328

@@ -341,7 +342,9 @@ Here, `d` specifies the dimension of each qudit. Default to `d=2` (qubit).
341342
"""
342343
function ghz_state(::Val{n}; d::Int = 2) where {n}
343344
nzind = collect((0:(d-1)) .* Int((d^n - 1) / (d - 1)) .+ 1)
344-
nzval = ones(ComplexF64, d) / sqrt(d)
345-
return QuantumObject(SparseVector(d^n, nzind, nzval), Ket(), ntuple(x -> d, Val(n)))
345+
nzval = fill(ComplexF64(1 / sqrt(d)), d)
346+
data = zeros(ComplexF64, d^n)
347+
@inbounds data[nzind] .= nzval
348+
return QuantumObject(data, Ket(), ntuple(x -> d, Val(n)))
346349
end
347350
ghz_state(n::Int; d::Int = 2) = ghz_state(Val(n), d = d)

0 commit comments

Comments
 (0)