DickeModel.TWA
For examples of usage, go to Examples for TWA.
Monte-Carlo systems
DickeModel.TWA.monte_carlo_integrate — Functionfunction 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
systemshould be an instance ofClassicalSystems.ClassicalSystem.mc_systemshould be an instance ofMonteCarloSystem.
Keyword arguments
tolerate_errorsindicates 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 than100errors occur consecutively, then then the procedure is stopped. Defaults totrue.show_progressis a boolean that toggles the progress bar. (Default istrue).maxNBatchis the maximum number of batch-sizes sent to each worker. Defaults toinf.kargsare redirected to [ClassicalSystems.integrate].
DickeModel.TWA.MonteCarloSystem — Typestruct MonteCarloSystemThis object may be passed to monte_carlo_integrate. Use mcs_for_averaging, mcs_for_variance, and mcs_for_survival_probability to generate them, and mcs_chain to join them together.
DickeModel.TWA.mcs_chain — Functionfunction 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.
function mcs_chain(mcs1::MonteCarloSystem, mcs2::MonteCarloSystem, ...)Functionally equivalent to mcs_chain([mcs1, mcs2, ...]).
Wigner distribution of coherent states
DickeModel.TWA.PhaseSpaceDistribution — Typestruct PhaseSpaceDistributionThis object represents a probability distribution in the phase space. Currently, the only implementation is through coherent_Wigner_SU2, coherent_Wigner_HW, and coherent_Wigner_HWxSU2
DickeModel.TWA.coherent_Wigner_HW — Functionfunction 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])
function coherent_Wigner_HW(u₀::AbstractVector{<:Real},j::Real=1,ħ::Real=1/j)Returns `coherentWignerHW(;q₀=u₀[1], p₀=u₀[2], j=j, ħ=ħ).
DickeModel.TWA.coherent_Wigner_SU2 — Functionfunction 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.
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, ħ=ħ).
DickeModel.TWA.coherent_Wigner_HWxSU2 — Functionfunction 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.
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, ħ=ħ).
Averages
DickeModel.TWA.mcs_for_averaging — Functionfunction 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
systemshould be an instance ofClassicalSystems.ClassicalSystem,
Keyword arguments
observablecan a function in the formf(x::Vector)orf(x1,x2 ,..., x_n)with n the dimension of the phase space determined bysystem.observablecan also be an expression determining the operations between thevarnamesdetermined by system. For example, ifsystemwere a instance ofClassicalDicke.DickeSystem, thenobservablecould be:(q+p^2 +Q),f(Q, q, P, p) = q + p^2 - Qorf(x) = x[2] + x[4]^2 - x[1], which are all equivalent. See also the submoduleWeyl, which produces expressions corresponding to the Weyl symbols of quantum observables. Note:observablecan also be an array of observables, in which case an array is returned for each time.tsshould be a sorted array of times.Nis the number of points to sample. The bigger the more accurate.distributionshould be an instance ofPhaseSpaceDistribution
DickeModel.TWA.average — Functionfunction average(system::ClassicalSystems.ClassicalSystem,
[[[same kargs as mcs_for_averaging]]],
kargs...)Calls mcs_for_averaging and then monte_carlo_integrate on the resulting MonteCarloSystem. Extra kargs are sent to the latter.
Variances
DickeModel.TWA.mcs_for_variance — Functionfunction 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.
DickeModel.TWA.variance — Functionfunction variance(system::ClassicalSystems.ClassicalSystem,
[[[same kargs as mcs_for_variance]]],
kargs...)Calls mcs_for_variance and then monte_carlo_integrate on the resulting MonteCarloSystem. Extra kargs are sent to the latter.
Survival probability
DickeModel.TWA.mcs_for_survival_probability — Functionfunction 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
systemshould be an instance ofClassicalSystems.ClassicalSystem,
Keyword arguments
Nis the number of points to sample. The bigger the more accurate.tsshould be a sorted array of times.distributionshould be an instance ofPhaseSpaceDistribution
DickeModel.TWA.survival_probability — Functionfunction survival_probability(system::ClassicalSystems.ClassicalSystem,
[[[same kargs as mcs_for_survival_probability]]],
kargs...)Calls mcs_for_survival_probability and then monte_carlo_integrate on the resulting MonteCarloSystem. Extra kargs are sent to the latter.
Matrix distributions
DickeModel.TWA.mcs_for_distributions — Functionfunction 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
xis:t, a matrix of dimensionslength(ts)$\times$length(ys)is produced. The coordinate(tᵢ,yⱼ)gives the probability of finding the observableybetweenys[j]andys[j+1]at the timets[i]. - If
yis:t, the same as above, changingxbyy. - If neither
xnoryare:t, andanimate = false, then a matrix of dimensionslength(xs)$\times$length(ys)is produced. The coordinate(xᵢ,yⱼ)gives the probability of finding the observablexbetweenxs[i]andxs[i+1]and the observableybetweenys[j]andys[j+1], averaging over all the times ints - If neither
xnoryaret, andanimate = truethen an array of matrices of dimensionslength(xs)$\times$length(ys)is produced, the same size asts. For thenth matrix in this array, the coordinate(xᵢ,yⱼ)gives the probability of finding the observablexbetweenxs[i]andxs[i+1]and the observableybetweenys[j]andys[j+1]at timets[n].
Arguments
systemshould be an instance ofClassicalSystems.ClassicalSystem,
Keyword arguments
xandycan be:t, or an observable as described in the arguments ofmcs_for_averaging.xsandysshould be given ifxandyare not:t, in wich case they should be range objects (e.g.0:0.1:1) containing the bins for the histogram.tsshould be a range object (e.g.0:0.1:1) of times.animateshould be aBool, which determines the behavior if neitherxnoryaret. Defaults totrueNis the number of points to sample. The bigger the more accurate.distributionshould be an instance ofPhaseSpaceDistribution
Note: If ts and y are both not passed, then y is set to :t and ts to 0:0.
DickeModel.TWA.calculate_distribution — Functionfunction calculate_distribution(system::ClassicalSystems.ClassicalSystem,
[[[same kargs as mcs_for_distributions]]],
kargs...)Calls mcs_for_distributions and then monte_carlo_integrate on the resulting MonteCarloSystem. Extra kargs are sent to the latter.
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.Jx — Methodfunction Jx(j::Real)Returns the Weyl symbol of the operator $\hat{J}_x$. (See p. 114 of Ref. [18])
DickeModel.TWA.Weyl.Jx² — Methodfunction Jx²(j::Real)Returns the Weyl symbol of the operator $\hat{J}_x^2$. (See bottom of p. 128 of Ref. [18])
DickeModel.TWA.Weyl.Jy — Methodfunction Jy(j::Real)Returns the Weyl symbol of the operator $\hat{J}_y$. (See p. 114 of Ref. [18])
DickeModel.TWA.Weyl.Jy² — Methodfunction Jy²(j::Real)Returns the Weyl symbol of the operator $\hat{J}_y^2$. (See bottom of p. 128 of Ref. [18])
DickeModel.TWA.Weyl.Jz — Methodfunction Jz(j::Real)Returns the Weyl symbol of the operator $\hat{J}_z$. (See p. 114 of Ref. [18])
DickeModel.TWA.Weyl.Jz² — Methodfunction Jz²(j::Real)Returns the Weyl symbol of the operator $\hat{J}_z^2$. (See bottom of p. 128 of Ref. [18])
DickeModel.TWA.Weyl.n — Methodfunction n(j::Real)Returns the Weyl symbol of the number operator $W(n̂)$, (p²+ q² - ħ)/2, where ħ = 1/j.
DickeModel.TWA.Weyl.n² — Methodfunction 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.