API Reference

SpectralDensities

SpectralDensities.AbstractSDType
abstract type AbstractSD

AbstractSD represents an abstract bath spectral density.

Any subtype of AbstractSD must at least define either sd or sdoverω for the new type.

source
SpectralDensities.DebyeSDType
struct DebyeSD <: AbstractSD

DebyeSD represents a Debye spectral density. It is characterized by an amplitude α representing the strength of the coupling and the cutoff frequency ωc. That is

\[J(\omega) = \frac{2\alpha}{\pi}\frac{\omega\omega_c^2}{\omega^2 + \omega_c^2}\]

Fields

  • α::Float64: The amplitude α, indicating the strength of the coupling.
  • ωc::Float64: The cutoff frequency.
source
SpectralDensities.DebyeSDMethod
DebyeSD(α::Real, ωc::Real)

Construct a Debye spectral density with the given amplitude α and cutoff frequency ωc.

Arguments

  • α::Real: The amplitude α, indicating the strength of the coupling.
  • ωc::Real: The cutoff frequency.

Returns

  • An instance of the DebyeSD struct representing the Debye spectral density.
source
SpectralDensities.ExponentialCutoffSDType
struct ExponentialCutoffSD{T <: AbstractSD} <: AbstractSD

ExponentialCutoffSD represents a spectral density with an exponential frequency cutoff. It is parameterized by the underlying spectral density J, and the frequency cutoff ωcutoff. That is

\[J_\mathrm{exp}(\omega) = J(\omega)e^{-\omega/\omega_c}\]

Fields

  • J::T: The underlying spectral density.
  • ωcutoff::Float64: The frequency cutoff.
source
SpectralDensities.ExponentialCutoffSDMethod
ExponentialCutoffSD(J::AbstractSD, ωcutoff)

Construct a spectral density with an exponential frequency cutoff using the underlying spectral density J and the given frequency cutoff ωcutoff.

Arguments

  • J::AbstractSD: The underlying spectral density.
  • ωcutoff: The frequency cutoff.

Returns

  • An instance of the ExponentialCutoffSD struct representing the spectral density with an exponential frequency cutoff.
source
SpectralDensities.GaussianCutoffSDType
struct GaussianCutoffSD{T <: AbstractSD} <: AbstractSD

GaussianCutoffSD represents a spectral density with a Gaussian frequency cutoff. It is parameterized by the underlying spectral density J, and the frequency cutoff ωcutoff. That is

\[J_\mathrm{gauss}(\omega) = J(\omega)e^{-\omega^2/\omega_c^2}\]

Fields

  • J::T: The underlying spectral density.
  • ωcutoff::Float64: The frequency cutoff.
source
SpectralDensities.GaussianCutoffSDMethod
GaussianCutoffSD(J::AbstractSD, ωcutoff)

Construct a spectral density with a Gaussian frequency cutoff using the underlying spectral density J and the given frequency cutoff ωcutoff.

Arguments

  • J::AbstractSD: The underlying spectral density.
  • ωcutoff: The frequency cutoff.

Returns

  • An instance of the GaussianCutoffSD struct representing the spectral density with a Gaussian frequency cutoff.
source
SpectralDensities.HardCutoffSDType
struct HardCutoffSD{T <: AbstractSD} <: AbstractSD

HardCutoffSD represents a spectral density with a hard frequency cutoff. It is parameterized by the underlying spectral density J, and the frequency cutoff ωcutoff. That is

\[J_\mathrm{hard}(\omega) = J(\omega)\Theta(\omega_c - \omega)\]

where $\Theta$ is the Heaviside theta function.

Fields

  • J::T: The underlying spectral density.
  • ωcutoff::Float64: The frequency cutoff.
source
SpectralDensities.HardCutoffSDMethod
HardCutoffSD(J::AbstractSD, ωcutoff)

Construct a spectral density with a hard frequency cutoff using the underlying spectral density J and the given frequency cutoff ωcutoff.

Arguments

  • J::AbstractSD: The underlying spectral density.
  • ωcutoff::Float64: The frequency cutoff.

Returns

  • An instance of the HardCutoffSD struct representing the spectral density with a hard frequency cutoff.
source
SpectralDensities.LorentzianSDType
struct LorentzianSD <: AbstractSD

LorentzianSD represents a Lorentzian spectral density. It is characterized by an amplitude α representing the strength of the coupling, the peak centre frequency ω0, and the width of the Lorentzian peak Γ. That is

\[J(\omega) = \frac{\alpha\Gamma}{\pi}\frac{\omega}{(\omega^2 - \omega_0^2)^2 + \omega^2\Gamma^2}\]

Fields

  • α::Float64: The amplitude α, indicating the strength of the coupling.
  • ω0::Float64: The centre frequency of the Lorentzian peak.
  • Γ::Float64: The width of the Lorentzian peak.
source
SpectralDensities.LorentzianSDMethod
LorentzianSD(α::Real, ω0::Real, Γ::Real)

Construct a Lorentzian spectral density with the given amplitude α, centre frequency ω0, and width Γ.

Arguments

  • α::Real: The amplitude α, indicating the strength of the coupling.
  • ω0::Real: The centre frequency of the Lorentzian peak.
  • Γ::Real: The width of the Lorentzian peak.

Returns

  • An instance of the LorentzianSD struct representing the Lorentzian spectral density.
source
SpectralDensities.OhmicSDType
struct OhmicSD <: AbstractSD

OhmicSD represents an Ohmic spectral density. It is characterized by an amplitude α representing the strength of the Ohmic coupling. That is

\[J(\omega) = \alpha\omega\]

Fields

  • α::Float64: The amplitude α, indicating the strength of the Ohmic coupling.
source
SpectralDensities.OhmicSDMethod
OhmicSD(α::Real)

Construct an Ohmic spectral density with amplitude α.

Arguments

  • α::Real: The amplitude α, indicating the strength of the Ohmic coupling.

Returns

  • An instance of the OhmicSD struct representing the Ohmic spectral density.
source
SpectralDensities.PolySDType
struct PolySD <: AbstractSD

PolySD represents a polynomial spectral density. It is characterized by an amplitude α representing the strength of the coupling and the polynomial degree n. That is

\[J(\omega) = \alpha\omega^n\]

Fields

  • α::Float64: The amplitude α, indicating the strength of the coupling.
  • n::Int: The polynomial degree.
source
SpectralDensities.PolySDMethod
PolySD(α::Real, n::Int)

Construct a polynomial spectral density with the given amplitude α and degree n.

Arguments

  • α::Real: The amplitude α, indicating the strength of the coupling.
  • n::Int: The polynomial degree.

Returns

  • An instance of the PolySD struct representing the polynomial spectral density.
source
SpectralDensities.correlationsMethod
correlations(J::AbstractSD, τ, β)

Calculate the correlation function for a spectral density J at a given time delay τ and inverse temperature β.

Arguments

  • J::AbstractSD: The spectral density.
  • τ: The time delay at which the correlation function is calculated.
  • β: The inverse temperature.

Returns

  • The correlation function for the spectral density J at the given time delay τ and inverse temperature β.
source
SpectralDensities.imag_memory_kernel_ftMethod
imag_memory_kernel_ft(J::AbstractSD, ω)

Calculate the imaginary part of the Fourier-transform of the memory kernel for a spectral density J at a given frequency ω.

Arguments

  • J::AbstractSD: The spectral density.
  • ω: The frequency at which the imaginary part of the Fourier-transform of the memory kernel is evaluated.

Returns

  • The imaginary part of the Fourier-transform of the memory kernel for the spectral density J at the given frequency ω.
source
SpectralDensities.memory_kernelMethod
memory_kernel(J::AbstractSD, τ)

Calculate the memory kernel for a spectral density J at a given time delay τ, that is

\[\mathcal{K}(\tau) = 2\Theta(\tau)\int_0^\infty J(\omega)\sin(\omega)\mathrm{d}\omega,\]

where $\Theta$ is the Heavisde theta function.

Arguments

  • J::AbstractSD: The spectral density.
  • τ: The time delay at which the memory kernel is evaluated.

Returns

  • The memory kernel for the spectral density J at the given time delay τ.
source
SpectralDensities.reorganisation_energyMethod
reorganisation_energy(J::AbstractSD)

Calculate the reorganization energy of a given spectral density J, i.e.

\[\int_0^\infty \frac{J(\omega)}{\omega} \mathrm{d}\omega \]

Arguments

  • J::AbstractSD: The spectral density.

Returns

  • The reorganization energy of the spectral density J.
source
SpectralDensities.sdMethod
sd(J::T, ω) where T <: AbstractSD

Evaluate the spectral density represented by J at a given frequency ω, i.e. J(ω).

Arguments

  • J::T: The spectral density.
  • ω: The frequency at which the spectral density is evaluated.

Returns

  • The spectral density J at the frequency ω.
source
SpectralDensities.sdoverωMethod
sdoverω(J::T, ω) where T <: AbstractSD

Evaluate the spectral density represented by J divided by a given frequency ω, i.e. J(ω)/ω.

Arguments

  • J::T: The spectral density.
  • ω: The frequency at which the spectral density is evaluated.

Returns

  • The spectral density J(ω) divided by ω.
source

WeakCoupling

SpectralDensities.WeakCoupling.cauchy_quadgkFunction
cauchy_quadgk(g, a, b, x0=0; kws...)

Computes the Cauchy principal value of the integral of a function g with singularity at x0 over the interval [a, b], that is

\[\mathrm{P.V.}\int_a^b\frac{g(x)}{x-x_0}\mathrm{d}x\]

Note that x0 must be contained in the interval [a, b]. The actual integration is performed by the quadgk method of the QuadGK.jl package and the keyword arguments kws are passed directly onto quadgk.

Note: If the function g contains additional integrable singularities, the user should manually split the integration interval around them, since currently there is no way of passing integration break points onto quadgk.

Arguments

  • g: The function to integrate.
  • a: The lower bound of the interval.
  • b: The upper bound of the interval.
  • x0: The location of the singularity. Default value is the origin.
  • kws...: Additional keyword arguments accepted by quadgk.

Returns

A tuple (I, E) containing the approximated integral I and an estimated upper bound on the absolute error E.

Throws

  • ArgumentError: If the interval [a, b] does not include the singularity x0.

Examples

julia> cauchy_quadgk(x -> 1/(x+2), -1.0, 1.0)
(-0.549306144334055, 9.969608472104596e-12)

julia> cauchy_quadgk(x -> x^2, 0.0, 2.0, 1.0)
(4.0, 2.220446049250313e-16)
source
SpectralDensities.WeakCoupling.hadamard_quadgkFunction
hadamard_quadgk(g, g′, a, b, x0=0; kws...)

Computes the Hadamard finite part of the integral of a function g with a singularity at x0 over the interval [a, b], that is

\[\mathcal{H}\int_a^b\frac{g(x)}{(x-x_0)^2}\mathrm{d}x\]

Note that x0 must be contained in the interval [a, b]. The actual integration is performed by the quadgk method of the QuadGK.jl package and the keyword arguments kws are passed directly onto quadgk.

Note: If the function g′ contains additional integrable singularities, the user should manually split the integration interval around them, since currently there is no way of passing integration break points onto quadgk.

Arguments

  • g: The function to integrate.
  • g′: The derivative of the function g.
  • a: The lower bound of the interval.
  • b: The upper bound of the interval.
  • x0: The location of the singularity. Defaults value is the origin.
  • kws...: Additional keyword arguments accepted by quadgk.

Returns

A tuple (I, E...) containing the approximated Hadamard finite part integral I and an additional error estimate E from the quadrature method.

Throws

  • ArgumentError: If the interval [a, b] does not include the singularity x0.

Examples

julia> hadamard_quadgk(x -> log(x+1), x -> 1/(x+1), 0.0, 2.0, 1.0)
(-1.6479184330021648, 9.969608472104596e-12)
source
SpectralDensities.WeakCoupling.weak_coupling_ΔMethod
weak_coupling_Δ(J::AbstractSD, ωB, β; ħ=one(ωB))

Calculate the weak-coupling coefficient Δ for the spectral density J, system Bohr frequency ωB, and inverse temperature β, defined as

\[\Delta_\beta(\omega_\mathrm{B}) = 2\omega_\mathrm{B}\int_0^\infty J(\omega)\frac{1}{\omega^2-\omega_\mathrm{B}^2} \coth\left(\frac{\beta\hbar\omega}{2}\right) \mathrm{d}\omega\]

See: J.D. Cresser, J. Anders, Phys. Rev. Lett. 127, 250601 (2021).

Arguments

  • J::AbstractSD: The spectral density.
  • ωB: The system Bohr frequency of interest.
  • β: The inverse temperature.
  • ħ: The value of the reduced Planck constant. Default is 1.

Returns

  • The weak-coupling coefficient Δ for the spectral density J, system Bohr frequency ωB, and inverse temperature β.
source
SpectralDensities.WeakCoupling.weak_coupling_ΔprimeMethod
weak_coupling_Δprime(J::AbstractSD, ωB, β; ħ=one(ωB))

Calculate the weak-coupling coefficient Δ′ for the spectral density J, system Bohr frequency ωB, and inverse temperature β, defined as

\[{\Delta'}_\beta(\omega_\mathrm{B}) = 2\int_0^\infty J(\omega)\frac{(\omega^2 + \omega_\mathrm{B}^2)}{(\omega^2-\omega_\mathrm{B}^2)^2} \coth\left(\frac{\beta\hbar\omega}{2}\right) \mathrm{d}\omega\]

See: J.D. Cresser, J. Anders, Phys. Rev. Lett. 127, 250601 (2021).

Note: The spectral density J must support automatic differentiation with ForwardDiff.jl.

Arguments

  • J::AbstractSD: The spectral density.
  • ωB: The system Bohr frequency of interest.
  • β: The inverse temperature.
  • ħ: The value of the reduced Planck constant. Default is 1.

Returns

  • The weak-coupling coefficient Δ′ for the spectral density J, system Bohr frequency ωB, and inverse temperature β.
source
SpectralDensities.WeakCoupling.weak_coupling_ΣMethod
weak_coupling_Σ(J::AbstractSD, ωB)

Calculate the weak-coupling coefficient Σ for the spectral density J and system Bohr frequency ωB, defined as

\[\Sigma(\omega_\mathrm{B}) = 2\int_0^\infty J(\omega)\frac{\omega}{\omega^2-\omega_\mathrm{B}^2} \mathrm{d}\omega\]

See: J.D. Cresser, J. Anders, Phys. Rev. Lett. 127, 250601 (2021).

Arguments

  • J::AbstractSD: The spectral denstiy.
  • ωB: The system Bohr frequency of interest.

Returns

  • The weak-coupling coefficient Σ for the spectral density J and system Bohr frequency ωB.
source
SpectralDensities.WeakCoupling.weak_coupling_ΣprimeMethod
weak_coupling_Σprime(J::AbstractSD, ωB)

Calculate the weak-coupling coefficient Σ′ for the spectral density J and system Bohr frequency ωB, defined as

\[\Sigma'(\omega_\mathrm{B}) = 4\omega_\mathrm{B}\int_0^\infty J(\omega)\frac{\omega}{(\omega^2-\omega_\mathrm{B}^2)^2} \mathrm{d}\omega\]

See: J.D. Cresser, J. Anders, Phys. Rev. Lett. 127, 250601 (2021).

Note: The spectral density J must support automatic differentiation with ForwardDiff.jl.

Arguments

  • J::AbstractSD: The spectral density.
  • ωB: The system Bohr frequency of interest.

Returns

  • The weak-coupling coefficient Σ′ for the spectral density J and system Bohr frequency ωB.
source