DickeModel.TWA

For examples of usage, go to Examples for TWA.

Monte-Carlo systems

DickeModel.TWA.monte_carlo_integrateFunction
function monte_carlo_integrate(system::ClassicalSystems.ClassicalSystem,
    mc_system::MonteCarloSystem;
    tolerate_errors::Bool = true,
    show_progress::Bool = true,
    maxNBatch::Real=Inf,kargs...)

This function is the backbone of this module. It performs a Monte Carlo integration-type procedure. The argument mc_system determines a distribution::PhaseSpaceDistribution, an integer N, and a list of times ts. This function calls an initializing function also determined by mc_system. Then, it samples N random points x from distribution, which are integrated using the Hamiltonian given by system. For each trajectory and time in ts, an operation is performed, which is determined again by mc_system. A final operation is then performed and the result is returned. This may seem a bit abstract, but it is a very flexible system. The applications are more concrete; for example, if you generate mc_system using mcs_for_averaging(... observable=f ...), then the initial operation is to generate an array of zeroes the same length as ts. Then, for each initial condition x, the result of f(x(t)) is added to each element of the array. Finally, the array is overall divided by N and then returned. This is exactly a Monte Carlo integration of the function f(x) over the classical evolution of the distribution.

This function uses all the available workers, but make sure to import this module in all of them.

Arguments

Keyword arguments

  • tolerate_errors indicates that some errors in the integration may be ignored. This is useful because sometimes a one-in-a-million numerical instability may arise, and you may want to ignore it. If more than 100 errors occur consecutively, then then the procedure is stopped. Defaults to true.
  • show_progress is a boolean that toggles the progress bar. (Default is true).
  • maxNBatch is the maximum number of batch-sizes sent to each worker. Defaults to inf.
  • kargs are redirected to [ClassicalSystems.integrate].
source
DickeModel.TWA.mcs_chainFunction
function mcs_chain(Fs::AbstractVector{MonteCarloSystem})

Generates a MonteCarloSystem by chaining together those in the array Fs = [monteCarloSystem1, monteCarloSystem2, ...]. The output of the system generated is an array that has the outputs of each system, in the same order as Fs.

If Fs may be any iterable. If it is not an AbstractVector, collect will be called on it.

source
function mcs_chain(mcs1::MonteCarloSystem, mcs2::MonteCarloSystem, ...)

Functionally equivalent to mcs_chain([mcs1, mcs2, ...]).

source

Wigner distribution of coherent states

DickeModel.TWA.coherent_Wigner_HWFunction
function coherent_Wigner_HW(;q₀::Real,p₀::Real,j::Real=1,ħ::Real=1/j)

Returns a PhaseSpaceDistribution corresponding to the two-dimensional Wigner function of a coherent state of the Heisenberg-Weyl algebra (i.e. a standard coherent state) centered at q₀,p₀. A value of ħ may be passed, or pass j to set ħ = 1/j. (See Eq. (B.1) of Ref. [19])

source
function coherent_Wigner_HW(u₀::AbstractVector{<:Real},j::Real=1,ħ::Real=1/j)

Returns `coherentWignerHW(;q₀=u₀[1], p₀=u₀[2], j=j, ħ=ħ).

source
DickeModel.TWA.coherent_Wigner_SU2Function
function coherent_Wigner_SU2(;Q₀::Real,P₀::Real,j::Real=1,ħ::Real=1/j)

Similar to coherent_Wigner_HW, but for a coherent state for the SU(2) algebra (a Bloch coherent state). The approximation given by Eq. (B.4) of Ref. [19] is used.

source
function coherent_Wigner_SU2(u₀::AbstractVector{<:Real},j::Real=1,ħ::Real=1/j)

Returns coherent_Wigner_SU2(Q₀=u₀[1], P₀=u₀[2], j=j, ħ=ħ).

source
DickeModel.TWA.coherent_Wigner_HWxSU2Function
function coherent_Wigner_HWxSU2(;Q₀::Real, q₀::Real,
                                P₀::Real, p₀::Real,
                                j::Real=1, ħ::Real=1/j)

Produces the Wigner function corresponding to the tensor product of coherent_Wigner_HW and coherent_Wigner_SU2.

source
function coherent_Wigner_HWxSU2(u₀::AbstractVector{<:Real},j::Real=1,ħ::Real=1/j)

Returns coherent_Wigner_HWxSU2(Q₀=u₀[1], q₀=u₀[2], P₀=u₀[3], p₀=u₀[4], j=j, ħ=ħ).

source

Averages

DickeModel.TWA.mcs_for_averagingFunction
function mcs_for_averaging(
    system::ClassicalSystems.ClassicalSystem;
    observable,
    ts::AbstractVector{<:Real}=[0.0],
    N::Integer,
    distribution::PhaseSpaceDistribution)

Generates a MonteCarloSystem that computes

\[ \frac{1}{N} \sum_{i=1}^N \text{observable}(\mathbf{x}_i(t))\]

for each $t$ in ts, where $\mathbf{x}_i$ will be sampled from distribution. The system produces an array the same size as ts, containing the result for each time.

Arguments

Keyword arguments

  • observable can a function in the form f(x::Vector) or f(x1,x2 ,..., x_n) with n the dimension of the phase space determined by system. observable can also be an expression determining the operations between the varnames determined by system. For example, if system were a instance of ClassicalDicke.DickeSystem, then observable could be :(q+p^2 +Q), f(Q, q, P, p) = q + p^2 - Q or f(x) = x[2] + x[4]^2 - x[1], which are all equivalent. See also the submodule Weyl, which produces expressions corresponding to the Weyl symbols of quantum observables. Note: observable can also be an array of observables, in which case an array is returned for each time.
  • ts should be a sorted array of times.
  • N is the number of points to sample. The bigger the more accurate.
  • distribution should be an instance of PhaseSpaceDistribution
source

Variances

DickeModel.TWA.mcs_for_varianceFunction
function mcs_for_variance(
    system::ClassicalSystems.ClassicalSystem;
    observable,
    ts::AbstractVector{<:Real}=[0.0],
    N::Integer,
    distribution::PhaseSpaceDistribution,
    return_average::Bool = false)

Generates a MonteCarloSystem that computes the variance

\[ \frac{1}{N} \sum_{i=1}^N \text{observable}^2(\mathbf{x}_i(t)) - \left (\frac{1}{N}\sum_{i=1}^N \text{observable}(\mathbf{x}_i(t))\right )^2\]

for each $t$ in ts, where $\mathbf{x}_i$ will be sampled from distribution. The system produces an array the same size as ts, containing the result for each time.

Arguments

  • See arguments and keyword arguments that mcs_for_averaging.
  • If return_average = false (default), only produces the variance, else it will produce a tuple, where the first element is the variance and the second the average.
source

Survival probability

DickeModel.TWA.mcs_for_survival_probabilityFunction
function mcs_for_survival_probability(
    system::ClassicalSystems.ClassicalSystem;
    N::Integer,
    ts::AbstractArray{<:Real},
    distribution::PhaseSpaceDistribution)

Generates a MonteCarloSystem that computes the survival probability through Eq. (C.7) of Ref. [19] (with $M=$ N, $w=$ distribution).

Arguments

Keyword arguments

  • N is the number of points to sample. The bigger the more accurate.
  • ts should be a sorted array of times.
  • distribution should be an instance of PhaseSpaceDistribution
source

Matrix distributions

DickeModel.TWA.mcs_for_distributionsFunction
function mcs_for_distributions(system::ClassicalSystems.ClassicalSystem;
    x,
    xs::Union{AbstractRange{<:Real},Nothing}=nothing,
    y=nothing,
    ts::Union{AbstractRange{<:Real},Nothing}=nothing,
    animate::Bool=false,
    ys::Union{AbstractRange{<:Real},Nothing}=nothing,
    N::Integer,
    distribution::PhaseSpaceDistribution)

Generates a MonteCarloSystem that produces a multidimensional histogram allowing to visualize the expected value of observables under PhaseSpaceDistribution distribution.

Several behaviors can be produced.

  • If x is :t, a matrix of dimensions length(ts) $\times$ length(ys) is produced. The coordinate (tᵢ,yⱼ) gives the probability of finding the observable y between ys[j] and ys[j+1] at the time ts[i].
  • If y is :t, the same as above, changing x by y.
  • If neither x nor y are :t, and animate = false, then a matrix of dimensions length(xs) $\times$ length(ys) is produced. The coordinate (xᵢ,yⱼ) gives the probability of finding the observable x between xs[i] and xs[i+1] and the observable y between ys[j] and ys[j+1], averaging over all the times in ts
  • If neither x nor y are t, and animate = true then an array of matrices of dimensions length(xs) $\times$ length(ys) is produced, the same size as ts. For the nth matrix in this array, the coordinate (xᵢ,yⱼ) gives the probability of finding the observable x between xs[i] and xs[i+1] and the observable y between ys[j] and ys[j+1] at time ts[n].

Arguments

Keyword arguments

  • x and y can be :t, or an observable as described in the arguments of mcs_for_averaging.
  • xs and ys should be given if x and y are not :t, in wich case they should be range objects (e.g. 0:0.1:1) containing the bins for the histogram.
  • ts should be a range object (e.g. 0:0.1:1) of times.
  • animate should be a Bool, which determines the behavior if neither x nor y are t. Defaults to true
  • N is the number of points to sample. The bigger the more accurate.
  • distribution should be an instance of PhaseSpaceDistribution

Note: If ts and y are both not passed, then y is set to :t and ts to 0:0.

source

Weyl symbols

The submodule TWA.Weyl generates classical phase-space observables that may be passed as the argument observable of average, calculate_distribution, etc. All the following functions return SymEngine objects, so they may be operated as if they were numbers.

DickeModel.TWA.Weyl.nMethod
function n(j::Real)

Returns the Weyl symbol of the number operator $W(n̂)$, (p²+ q² - ħ)/2, where ħ = 1/j.

source
DickeModel.TWA.Weyl.n²Method
function n²(j::Real)

Returns the Weyl symbol of the number operator squared $W(n̂^2) = W(n̂)^2 - \hbar^2/4$, where $\hbar = 1/j$. (Note that the extra term $\hbar^2/4$ appears due to the non-commutativity of $q̂$ and $p̂$. See Ref. [16] for details on how to compute these expressions.

source