DickeModel.DickeBCE

For examples of usage, go to Examples for DickeBCE.

DickeModel.DickeBCE.QuantumDickeSystemType
mutable struct QuantumDickeSystem <: ClassicalDicke.DickeSystem

This object represents the quantum Dicke model. It stores the parameters of the system, and it may be passed to multiple functions in this module. To generate it, use

function QuantumDickeSystem(classical_system::ClassicalDicke.ClassicalDickeSystem;
    j::Real,
    Nmax::Union{Integer,Nothing}=nothing)

or DickeBCE.QuantumDickeSystem(;ω₀,ω,γ,j,Nmax).

Because this object subtypes ClassicalDicke.DickeSystem, it may be passed to any function that requires a DickeSystem or a ClassicalSystem, and the underlying classical_system will be used.

Arguments

  • classical_system should be generated with ClassicalDicke.ClassicalDickeSystem.
  • j is the value of $j$. It must be a positive half-integer.
  • Nmax - 1 is the maximum excitation of the modified bosonic sector in the efficient coherent basis (see [4], [5]). Nmax must be a positive integer. It may be omitted. In this case, you may use functions that do not require Nmax. Moreover, if you call diagonalization with a system that has Nmax = nothing, and there is a diagonalization saved in disk, diagonalization(system) will load the eigenstates with the largest Nmax and will set the value of Nmax in system.
source
DickeModel.DickeBCE.dimensionFunction
function dimension(system::QuantumDickeSystem)

Returns the integer $(2j + 1)\times N_\text{max}$, which gives the dimension of the Hilbert space.

source

Diagonalization

DickeModel.DickeBCE.diagonalizationFunction
function diagonalization(system::QuantumDickeSystem;
        load_cache::Bool = true,
        save_cache::Bool = true,
        cache_folder::AbstractString = joinpath(homedir(),"dicke_diagonalizations"),
        maxϵ::Real = 5.0,
        onlyload::Union{AbstractVector{<:Integer},Nothing} = nothing,
        only_eigenenergies::Bool = false,
        verbose::Bool = true,
        converged_tolerance::Real = 1e-3)

Diagonalizes the Dicke Hamiltonian up to a maximum energy maxϵ. The resulting eigenstates are guaranteed to be converged, with a tolerance determined by converged_tolerance. Numerical degeneracies are also corrected, to ensure that the eigenstates have parity +1 or -1. If only_eigenenergies is false (default), a tuple (eigenenergies,eigenstates) is returned, where eigenenergies is real vector containing the eigenenergies and eigenstates is a real matrix that contains the eigenstates as columns.

If load_cache = true it will try to load saved diagonalizations from cache_folder. In the case that system has Nmax = nothing, the diagonalization with the largest Nmax will be loaded, and Nmax will be set in system.

Arguments

Keyword arguments

  • load_cache is a boolean indicating whether to try to load from cache folder. Defaults to true.
  • save_cache determines if the results of a computed diagonalization are saved to the cache folder. Defaults to true.
  • cache_folder is the cache folder where diagonalizations are saved. Default is %HOME%/dicke_diagonalizations
  • maxϵ is the maximum energy up to which we diagonalize. Keep this number higher than the maximum converged regime you want.
  • onlyload may be a vector of integers, indicating the indices of the eigenstates to load, or the default value, nothing, indicates that all eigenstates should be loaded.
  • only_eigenenergies should be true if you only to load the eigenenergies. Defaults to false (return eigenstates and eigenenergies).
  • verbose is true (default) if you want to see info messages.
  • converged_tolerance determines how strict we are in saying that an eigenstate is converged (see [5]). The default value 1e-3 is usually the best.
source

Evolution

DickeModel.DickeBCE.evolveFunction
function evolve(t::Union{Real,AbstractArray{<:Real}}
    state::AbstractVector{<:Number};
    eigenstates::AbstractMatrix{<:Number},
    eigenenergies::AbstractVector{<:Real},
    normwarning::Real=0.99)

Computes the evolution $e^{-i \hat{H} t}|$ state $\rangle$ under the Dicke Hamiltonian $\hat{H}$. The result is returned in the BCE.

Arguments

  • t may be a single time, or a vector of times. In the first case, a vector will be returned, and in the second, a matrix with each column corresponding to each time.
  • state is a vector representing the state in the BCE.

Keyword arguments

  • eigenstates and eigenenergies are the ones returned by diagonalization.
  • normwarning is a tolerance (defaults to 0.99). If the norm of the state in the eigenbasis is below this number, a warning will be thrown. This usually happens if Nmax is too small for the energy regime you are working on.
source
DickeModel.DickeBCE.survival_probabilityFunction
function survival_probability(t::Union{Real,AbstractArray{<:Real}};
    state::AbstractVector{<:Number},
    eigenstates::AbstractMatrix{<:Number},
    eigenenergies::AbstractVector{<:Real},
    normwarning::Real=0.99)

Computes the survival probability $| \langle$ state $| e^{-i \hat{H} t} |$ state $\rangle | ^2$ at the times given.

Arguments

  • t may be a single time, or a vector of times. In the first case, a number will be returned, and in the second, a vector.

Keyword arguments

  • state is a vector representing the state in the BCE.
  • eigenstates and eigenenergies are the ones returned by diagonalization.
  • normwarning is a tolerance (defaults to 0.99). If the norm of the state in the eigenbasis is below this number, a warning will be thrown. This usually happens if Nmax is too small for the energy regime you are working on.
source

Operators

Participation Ratio

DickeModel.DickeBCE.participation_ratioFunction
function participation_ratio(state::AbstractVector{<:Number};
    eigenstates::AbstractMatrix{<:Number},
    eigenenergies::Union{AbstractVector{<:Real},Nothing}=nothing,
    count_degeneracies::Bool=eigenenergies!=nothing,
    degentol::Real=1e-5)

Returns the Participation Ratio (PR) of state in the eigenbasis. If count_degeneracies is true, then Eq. (7) of [19] is used. If count_degeneracies is false, then the inverse of Eq. (19) of [19] is used.

Arguments

  • state should be a complex vector representing the state in the efficient coherent basis.

Keyword arguments

  • eigenstates should be a matrix containing the eigenstates.
  • eigenenergies should be passed if count_degeneracies is true. It is a list containing the eigenenergies.
  • count_degeneracies – whether to modify the PR definition to account for degeneracies (see above). Default is false if eigenenergies is not passed, else it is true.
  • degentol minimum energy separation below which two eigenstates are considered degenerate. Default is 1e-5.
source
DickeModel.DickeBCE.factor_R_of_coherent_stateFunction
function factor_R_of_coherent_state(system::QuantumDickeSystem,
        x::AbstractVector{<:Real};
        eigenstates::AbstractMatrix{<:Number},
        eigenenergies::AbstractVector{<:Real},
        state::AbstractVector{<:Number} = coherent_state(system,x))

Computes $R$ for a coherent state, as defined in Eq. (30) of [19].

Arguments

Keyword arguments

  • eigenstates should be a matrix containing the eigenstates as columns.
  • eigenenergies is a list containing the eigenenergies.
  • state is a complex vector representing the coherent state $\left | \mathbf{x} \right \rangle$. If it is not passed, it is computed using DickeBCE.coherent_state
source

Coherent states

DickeModel.DickeBCE.coherent_state!Function
function coherent_state!(system::QuantumDickeSystem,
            x::AbstractArray{<:Real,1},
            data::AbstractVector{<:Number};
            normwarning=0.99,
            add_state_to_data::Bool=false,
            extra_phase::Complex=1.0+0im,
            chop::Real=1e-6)

This function computes the coefficients of a coherent state centered at x in the BCE, and stores the result in data. (See Ref. [19])

Arguments

Keyword arguments

  • normwarning is a tolerance (defaults to 0.99). If the norm of the resulting coherent state is below this number, a warning will be thrown. This usually happens if Nmax is too small for the energy regime you are working on.
  • add_state_to_data is a boolean. If it is true, the coefficients will be added to data. If it is false (default), this function will override any preexisting values in data. This is useful if you want to add many coherent states together without having to allocate that much memory.
  • extra_phase is a complex number that multiplies the resulting state overall. Defaults to 1.
  • chop is a numerical tolerance between 0 and 1. If chop=0, then all the coefficients all computed. However, if ε = 1 - chop > 0.0, then some coefficients at the tail of the distribution (whose total squared norm does not exceed ε) will be treated as 0.0. This allows to significantly reduce computation time, but introduces a numerical error of order ε. The default is 1e-6. See this example and Ref. [11].
source
DickeModel.DickeBCE.coherent_stateFunction
function coherent_state(system::QuantumDickeSystem,
    x::AbstractVector{<:Real};
    normwarning::Real=0.99,
    extra_phase::Complex=1.0+0im,
    chop::Real=0.0)

Calls coherent_state! passing data as a new vector and returns the result.

source

Overlaps

DickeModel.DickeBCE.coherent_overlapFunction
function coherent_overlap(system::QuantumDickeSystem,
    x::AbstractVector{<:Real},
    state::AbstractVector{<:Number};
    chop::Real=1e-6,
    normwarning::Real=0.99,
    datacache::Union{AbstractArray{Complex{Float64},1},Nothing}=nothing)

Returns the overlap $\langle \mathbf{x} |$ state $\rangle$, where $\left | \mathbf{x} \right \rangle$ is a coherent state centered at x.

Arguments

Keyword arguments

  • chop is a numerical tolerance between 0 and 1. If chop = 0, then all the products of coefficients are computed. However, if ε = 1 - chop > 0.0, then some coefficients of the coherent state at the tail of the distribution (whose total squared norm does not exceed ε) will be treated as 0.0. This allows to significantly reduce computation time, but introduces a numerical error of order ε. The default is 1e-6. See this example and Ref. [11].
  • normwarning is a tolerance (defaults to 0.99). If the norm of the coherent state is below this number, a warning will be thrown. This usually happens if Nmax is too small for the energy regime you are working on.
  • datacache in an array where to store the result, or nothing for the result to be returned.
source
function coherent_overlap(system::QuantumDickeSystem,
    x::AbstractVector{<:Real},
    states::AbstractMatrix{<:Number};
    chop::Real=1e-6,
    normwarning::Real=0.99,
    datacache::Union{AbstractVector{<:Complex},Nothing}=nothing)

Same as coherent_overlap(..., state,...), but instead of one state, it allows states to be a matrix with multiple states as columns. It is functionally equivalent to

[coherent_overlap(..., states[:,k], ...) for k in 1:size(states)[2]]

but it is (much) faster. The result is stored in datacache if provided.

source
DickeModel.DickeBCE.husimi_of_coherentFunction
function husimi_of_coherent(system::QuantumDickeSystem,
    x::AbstractArray{<:Real,1},
    y::AbstractArray{<:Real,1})

Returns $\left | \left \langle \mathbf{x}\middle | \mathbf{y}\right \rangle \right | ^2$. Equation (3.14b) of Ref. [2] is used for the overlap of the Bloch coherent states.

Arguments

source

Semiclassical approximations

DickeModel.DickeBCE.density_of_statesFunction
function density_of_states(system::QuantumDickeSystem, ϵ::Real)

Returns the semiclassical density of states (DoS) $ν(ϵ)$, in units of $1/ϵ$, as given by Eq. (A1) of Ref. [15].

This is computed by multiplying the volume of the classical energy shell at ϵ by $(2\pi \hbar_\text{eff})^2=(2 \pi / j)^2$ (see ClassicalDicke.energy_shell_volume).

Arguments

See Plotting the semiclassical density of states for an example.

source
DickeModel.DickeBCE.energy_width_of_coherent_stateFunction
function energy_width_of_coherent_state(system::QuantumDickeSystem,
                                        x::AbstractVector{<:Real})

Returns the energy width $\sigma$ of the coherent state $\left | \mathbf{x}\right \rangle$, in units of $\epsilon$. This quantity is given by $\sigma_D/j$ with $\sigma_D$ as in App. A of Ref. [7].

Arguments

  • system should be an instance of DickeBCE.QuantumDickeSystem.
  • x is the coordinate $\mathbf{x}$ of the coherent state in the format [Q,q,P,p].
source

Random states

DickeModel.DickeBCE.random_cₖ_generatorFunction
function random_cₖ_generator(system::QuantumDickeSystem;
                            kargs...)

Generates coefficients to build random states, which are given by

\[ c_k(\epsilon_k) \sim e^{i \theta_k} \sqrt{\frac{r_k\rho(\epsilon_k)}{\nu(\epsilon_k)}}.\]

(See Refs. [8], [19], and [14]). It returns a function cₖ(ϵ;number=1,cache=nothing), which produces n=number coefficients, storing them in cache (an n-vector), if given.

The parameters of the equation above are determined as follows:

Required keyword arguments

  • To determine $\rho(\epsilon)$ above, pass either:
    • σ and ϵ, both real, in which case $\rho$ will be taken as a normal distribution centered at ϵ with standard deviation σ,
    • OR envelope as any of Distributions.UnivariateDistribution.
  • To determine $r_k$, pass either:
    • ensemble = :GUE, which takes $r_k$ from an exponential distribution Exponential(0.91) or ensemble = :GOE, which takes $r_k$ from a $\chi^2$ distribution with one degree of freedom,
    • OR rₖ_distribution as any of Distributions.UnivariateDistribution.
  • To determine $\theta_k$, pass either:
    • phases = :real, which takes $\theta_k$ to be random from $\{ -\pi, \pi \}$ or phases = :complex, which takes $\theta_k$ to be uniform in $[0, 2\pi)$,
    • OR ensemble = :GOE, which automatically sets phases = :real, or ensemble = :GUE, which automatically sets phases = :complex.

Optional keyword arguments

  • divide_by_DoS is a boolean. If true (default), leaves the density of states, $\nu(\epsilon_k)$ in the equation above. A value of false removes it.
source
DickeModel.DickeBCE.random_state!Function
function random_state!(system::QuantumDickeSystem,
    data::AbstractVecOrMat{<:Number};
    kargs...)

Generates as many random states as columns in data saving the result to data. The general form of the random states in the eigenbasis is $\left | R \right \rangle = \sum_k c_k \left | E_k \right \rangle$, where the $c_k$ are computed using random_cₖ_generator. You may choose to select only positive or negative parity eigenstates (see below).

Arguments

  • system should be an instance of DickeBCE.QuantumDickeSystem.
  • data should be a vector or matrix of numbers, where the result is stored. It must have as many columns as random numbers you want and dimension(system) rows.

Required keyword arguments

Optional keyword arguments

  • All the optional keyword arguments of random_cₖ_generator may be passed.
  • tol may be a real number which determines the tolerance for the state convergence. A value tol=0 builds the state with all of the eigenstates, and a positive number cuts the tails of the energy envelope up to that number. Higher values are faster, but decrease precision (Defaults to 1e-6)
  • parity may be +, -, or nothing (default). If it is not nothing, the random will be composed of only the eigenstates with this parity.
  • parities should be passed if parity is not nothing. It should be a vector of -1s and 1s containing the parities of all of the eigenstates. If it is not passed, it will be computed with eigenstate_parities, but this is slow. If you plan to use this function repeatedly, you should precompute parities by calling eigenstate_parities yourself and then pass it to parities to avoid repeated calls to eigenstate_parities.
source
DickeModel.DickeBCE.random_stateFunction
function random_state(system::QuantumDickeSystem,
                    n::Int=1;
                    kargs...)

Generates n random states, returning them as columns in a matrix (or a vector if n = 1, which is the default). kargs are redirected to random_state!: some additional keyword arguments are required: see random_state!.

source

Wigner functions

DickeModel.DickeBCE.WignerFunction
function Wigner(system::QuantumDickeSystem,
    state::AbstractArray{<:Number},
    points::Vector{<:AbstractVector{<:Real}})

Evaluates the Wigner function of state (which is a vector in the BCE) in all points = [[Q1,q1,P1,p1], [Q2,q2,P2,p2], ...], returning an array with the results.

Note: This function has not been thoroughly tested. It is based on Ref. [10].

source
DickeModel.DickeBCE.WignerProjQPFunction
function WignerProjQP(system::QuantumDickeSystem,
    state::AbstractArray{<:Number},
    ptsqp::Vector{<:AbstractVector{<:Real}})

Evaluates Wigner function, integrated over the atomic variables $Q,P$, of state in all ptsqp = [[q1,p1], [q2,p2], ...], returning an array with the results.

Note: This function has not been thoroughly tested. It is based on Ref. [10].

source
DickeModel.DickeBCE.WignerProjqpFunction
function WignerProjqp(system::QuantumDickeSystem,
    states::AbstractArray{<:AbstractArray{<:Number}},
    ptsQP::Vector{<:AbstractVector{<:Real}};
    show_progress::Bool=true)

Evaluates Wigner function, integrated over the bosonic variables $q,p$, of each of the states in states (which is a vector of complex vectors) in all ptsQP = [[Q1,P1], [Q2,P2], ...], returning an array with the results. See this example. If show_progress = true it shows a progress bar.

Note: This function has not been thoroughly tested. It is based on Ref. [10].

source