Dynamics
SpiDy.diffeqsolver
— Functionfunction diffeqsolver(s0, tspan, J::LorentzianSD, Jshared::LorentzianSD, bfields, bfieldshared, matrix::Coupling; JH=zero(I), S0=1/2, Bext=[0, 0, 1], saveat=[], save_fields=false, projection=true, alg=Tsit5(), atol=1e-3, rtol=1e-3)
Solves the dynamics of a system of ineracting spins under the influence of both local (unique to each spin) and global (shared by all spins) stochastic noise from the environment.
Arguments
s0
: Array of length3N
specifying the initial conditions of theN
spins. The order the initial consitions is first theSx,Sy,Sz
for the first spin, then for the second, and so on.tspan
: The time span to solve the equations over, specified as a tuple(tstart, tend)
.Jlist::Vector{LorentzianSD}
: The list of spectral densities of the environment(s).bfield
: An array of tuples of functionsArray{Tuple{Function, Function, Function}}
representing the time series of the noise for each environment.bcoupling::Vector{Vector{T}} T <: Real
: A vector (of length the number of baths) of vectors (of length the number of spins), specifying if each bath couples to each spin.matrix::Vector{TT} TT <: Coupling
: A vector ofCoupling
structs specifying the spin-environment coupling matrix for each environment.JH=zero(I)
: (Optional) The spin-spin coupling matrix. Default is zero matrix (i.e. non-interacting spins).S0=1/2
: (Optional) The spin quantum number. Default is 1/2.Bext=[0, 0, 1]
: (Optional) The external magnetic field vector. Default is[0, 0, 1]
(normalised length pointing in thez
direction).saveat=[]
: (Optional) An array of time points where the solution should be saved. Default is empty, which saves the solution at the time steps chosen by the integration algorithm.save_fields=false
: (Optional) If true, also return the auxiliary fields encoding the environment memory.projection=false
: (Optional) Specifies whether to project the spin vectors onto the unit sphere at each time step, hence forcing the numerical conservation of the spin length. Default isfalse
.alg=Vern6()
: (Optional) The differential equation solver algorithm. Default isVern6()
. See theDifferentialEquations.jl
docs for choices.atol=1e-6
: (Optional) The absolute tolerance for the solver. Default is1e-6
.rtol=1e-6
: (Optional) The relative tolerance for the solver. Default is1e-6
.
Note: The LorentzianSD
type is provided by the SpectralDensities.jl package.
Note: Additional keyword arguments will be passed on to the ODE solver (see the DifferentialEquations.jl
docs)
Returns
An ODESolution
struct from DifferentialEquations.jl
containing the solution of the equations of motion.
function diffeqsolver(s0, tspan, J::LorentzianSD, Jshared::LorentzianSD, bfields, bfieldshared, matrix::Coupling; JH=zero(I), S0=1/2, Bext=[0, 0, 1], saveat=[], save_fields=false, projection=true, alg=Tsit5(), atol=1e-3, rtol=1e-3)
Solves the dynamics of a system of ineracting spins under the influence of either local (unique to each spin) or global (shared by all spins) stochastic noise from the environment.
Arguments
s0
: Array of length3N
specifying the initial conditions of theN
spins. The order the initial consitions is first theSx,Sy,Sz
for the first spin, then for the second, and so on.tspan
: The time span to solve the equations over, specified as a tuple(tstart, tend)
.J::LorentzianSD
: The spectral density of the noise acting on the spins (either local or shared depending on the value ofbfields
).bfields
: For local baths, an array of tuples of functionsArray{Tuple{Function, Function, Function}}
representing the time series of the local stochastic field for each spin. For a global bath, a tuple of functionsTuple{Function, Function, Function}
representing the time series of the global stochastic field shared by all the spins.matrix::Coupling
: The spin-environment coupling matrix.JH=zero(I)
: (Optional) The spin-spin coupling matrix. Default is zero matrix (i.e. non-interacting spins).S0=1/2
: (Optional) The spin quantum number. Default is 1/2.Bext=[0, 0, 1]
: (Optional) The external magnetic field vector. Default is[0, 0, 1]
(normalised length pointing in thez
direction).saveat=[]
: (Optional) An array of time points where the solution should be saved. Default is empty, which saves the solution at the time steps chosen by the integration algorithm.save_fields=false
: (Optional) If true, also return the auxiliary fields encoding the environment memory.projection=false
: (Optional) Specifies whether to project the spin vectors onto the unit sphere at each time step, hence forcing the numerical conservation of the spin length. Default isfalse
.alg=Vern6()
: (Optional) The differential equation solver algorithm. Default isVern6()
. See theDifferentialEquations.jl
docs for choices.atol=1e-6
: (Optional) The absolute tolerance for the solver. Default is1e-6
.rtol=1e-6
: (Optional) The relative tolerance for the solver. Default is1e-6
.
Note: The LorentzianSD
type is provided by the SpectralDensities.jl package.
Note: Additional keyword arguments will be passed on to the ODE solver (see the DifferentialEquations.jl
docs)
Returns
An ODESolution
struct from DifferentialEquations.jl
containing the solution of the equations of motion.
function diffeqsolver(x0, p0, tspan, J::LorentzianSD, bfield, matrix::Coupling; JH=zero(I), Ω=1.0, counter_term=true, saveat=[], save_fields=false, alg=Tsit5(), atol=1e-3, rtol=1e-3)
Solves the dynamics of a system of ineracting oscillators under the influence of either local (unique to each oscillator) or global (shared by all oscillators) stochastic noise from the environment.
Arguments
x0
: Array of lengthN
specifying the initial position of theN
oscillators.p0
: Array of lengthN
specifying the initial momentum of theN
oscillators.tspan
: The time span to solve the equations over, specified as a tuple(tstart, tend)
.J::LorentzianSD
: The spectral density of the noise acting globally on the harmonic oscillators.bfield
: A tuple of functionsTuple{Function, Function, Function}
representing the time series of the global stochastic field shared by all the harmonic oscillators.matrix::Coupling
: The harmonic oscillators-environment coupling matrix.JH=zero(I)
: (Optional) The oscillator-oscillator coupling matrix. Default is zero matrix (i.e. non-interacting oscillators).Ω=1.0
: (Optional) The natural angular frequency of the harmonic oscillators (currently the same for all). Default is 1.counter_term=true
: (Optional) Whether to include the counter-term or not. Default is true.saveat=[]
: (Optional) An array of time points where the solution should be saved. Default is empty, which saves the solution at the time steps chosen by the integration algorithm.save_fields=false
: (Optional) If true, also return the auxiliary fields encoding the environment memory.alg=Vern6()
: (Optional) The differential equation solver algorithm. Default isVern6()
. See theDifferentialEquations.jl
docs for choices.atol=1e-6
: (Optional) The absolute tolerance for the solver. Default is1e-6
.rtol=1e-6
: (Optional) The relative tolerance for the solver. Default is1e-6
.
Note: The LorentzianSD
type is provided by the SpectralDensities.jl package.
Note: Additional keyword arguments will be passed on to the ODE solver (see the DifferentialEquations.jl
docs)
Returns
An ODESolution
struct from DifferentialEquations.jl
containing the solution of the equations of motion.