Welcome!

KadanoffBaym.jl is the first fully adaptive solver for Kadanoff-Baym equations written in Julia.

Tip

To learn more about the solver and Kadanoff-Baym equations, have a look into our accompanying paper.

Installation

To install, use Julia's built-in package manager

julia> ] add KadanoffBaym

Scalability

For now, KadanoffBaym.jl is restricted to run on a single machine, for which the maximum number of threads available will be used. You can set this number by running Julia with the thread flag

julia -t auto

Examples

To learn how to work with KadanoffBaym.jl, there are two options:

  • The examples folder of our repository, which contains notebooks for all of the systems studied in our paper.

  • The examples section of this documentation. If you are interested in quantum dynamics, we recommend you start with the tight-binding model. More advanced users can jump directly to Fermi-Hubbard model part I about the second Born approximation. Part II shows how to solve the more involved $T$-matrix approximation.

Note

KadanoffBaym.jl can also be used to simulate stochastic processes. An introduction to this topic is given here.

Library

KadanoffBaym.jl was designed to be lean and simple and hence only exports a handful of functions, namely GreenFunction (together with two possible time symmetries, Symmetrical and SkewHermitian) and the integrator kbsolve!. Besides these, wigner_transform can be used to analyze data in a Wigner(-Ville) transformed basis and TimeOrderedGreenFunction to quickly compute the Keldysh components (greater, lesser, retarded and advanced) of time convolutions via the Langreth rules.

Note

You need to import an FFT library – e.g., FFTW – to use wigner_transform.

Index

Solver

KadanoffBaym.kbsolve!Function
kbsolve!(fv!::Function, fd!::Function, u0::Vector{<:AbstractGreenFunction}, (t0, tmax)::Tuple{Union{Real, Vector{<:Real}}, Real})

Solves the 2-time Voltera integro-differential equation

$du(t_1,t_2) / dt_1 = f_v(t_1,t_2) = v[u,t_1,t_2] + ∫_{t0}^{t1} dτ K_1^v[u,t_1,t_2,τ] + ∫_{t0}^{t2} dτ K_2^v[u,t_1,t_2,τ]$

$du(t_1,t_2) / dt_2 = f_h(t_1,t_2) = h[u,t_1,t_2] + ∫_{t0}^{t1} dτ K_1^h[u,t_1,t_2,τ] + ∫_{t0}^{t2} dτ K_2^h[u,t_1,t_2,τ]$

for 2-point functions u0 from t0 to tmax.

Parameters

  • fv!(out, ts, w1, w2, t1, t2): The right-hand side of $du(t_1,t_2)/dt_1$ on the time-grid (ts x ts). The weights w1 and w2 can be used to integrate the Volterra kernels K1v and K2v as sum_i w1_i K1v_i and sum_i w2_i K2v_i, respectively. The output is saved in-place in out, which has the same shape as u0
  • fd!(out, ts, w1, w2, t1, t2): The right-hand side of $(du(t_1,t_2)/dt_1 + du(t_1,t_2)/dt_2)|_{t_2 → t_1}$
  • u0::Vector{<:GreenFunction}: List of 2-point functions to be integrated
  • (t0, tmax): A tuple with the initial time(s) t0 – can be a vector of past times – and final time tmax

Optional keyword parameters

  • f1!(out, ts, w1, t1): The right-hand-side of $dv(t_1)/dt_1$. The weight w1 can be used to integrate the Volterra kernel and the output is saved in-place in out, which has the same shape as v0
  • v0::Vector{<:GreenFunction}: List of 1-point functions to be integrated
  • callback(ts, w1, w2, t1, t2): A function that gets called everytime the 2-point function at indices (t1, t2) is updated. Can be used to update functions which are not being integrated, such as self-energies
  • stop(ts): A function that gets called at every time-step that stops the integration when it evaluates to true
  • atol::Real: Absolute tolerance (components with magnitude lower than atol do not guarantee number of local correct digits)
  • rtol::Real: Relative tolerance (roughly the local number of correct digits)
  • dtini::Real: Initial step-size
  • dtmax::Real: Maximal step-size
  • qmax::Real: Maximum step-size factor when adjusting the time-step
  • qmin::Real: Minimum step-size factor when adjusting the time-step
  • γ::Real: Safety factor for the calculated time-step such that it is accepted with a higher probability
  • kmax::Integer: Maximum order of the adaptive Adams method
  • kmax_vie::Integer: Maximum order of interpolant of the Volterra integrals Heuristically, it seems that having too high of a kmax_vie can result in numerical instabilities

Notes

  • Due to high memory and computation costs, kbsolve! mutates the initial condition u0 and only works with in-place rhs functions, unlike standard ODE solvers.
  • The Kadanoff-Baym timestepper is a 2-time generalization of the variable Adams method presented in E. Hairer, S. Norsett and G. Wanner, Solving Ordinary Differential Equations I: Non- stiff Problems, vol. 8, Springer-Verlag Berlin Heidelberg, ISBN 978-3-540-56670-0, doi:10.1007/978-3-540-78862-1 (1993).
source

Green Functions

KadanoffBaym.GreenFunctionType
GreenFunction(g::AbstractArray, s::AbstractSymmetry)

A container interface for g with array indexing respecting some symmetry rule s. Because of that, g must be square in its last 2 dimensions, which can be resized with resize!.

The array g is not restricted to being contiguous. For example, g can have Matrix{T}, Array{T,4}, Matrix{SparseMatrixCSC{T}}, etc as its type.

Notes

The GreenFunction does not own g. Proper care must be taken when using multiple GreenFunctions since using the same array will result in unexpected behaviour

julia> data = zeros(2,2)
julia> g1 = GreenFunction(data, Symmetrical)
julia> g2 = GreenFunction(data, Symmetrical)
julia> g1[1,1] = 3
julia> @show g2[1,1]
julia> g1.data === g2.data # they share the same data

Indexing with less indices than the dimension of g results in a "take-all-to-the-left" indexing

julia> gf[i,j] == gf[:,:,...,:,i,j]
julia> gf[i,j,k] == gf[:,:,...,:,i,j,k]

In order to ensure a correct behaviour of the symmetries, setindex! is only defined at the level of time coordinates

julia> gf[i, j] = a # is equivalent to gf[..., i, j] = a[...]

Custom symmetries can be implemented via multiple dispatch

julia> struct MySymmetry <: KadanoffBaym.AbstractSymmetry end
julia> @inline KadanoffBaym.symmetry(::Type{MySymmetry}) = conj

Examples

GreenFunction simply takes some data g and embeds the symmetry s in its indexing

julia> time_dim = 3
julia> spin_dim = 2
julia> data = zeros(spin_dim, spin_dim, time_dim, time_dim)
julia> gf = GreenFunction(data, Symmetrical)
julia> gf[2,1] = rand(spin_dim, spin_dim)
julia> @show gf[1,2]
julia> @show KadanoffBaym.symmetry(Symmetrical)(gf[2,1])
source

Wigner Transformation

KadanoffBaym.wigner_transformFunction
wigner_transform(x::AbstractMatrix; ts=1:size(x,1), fourier=true)

Wigner-Ville transformation

$x_W(ω, T) = i ∫dt x(T + t/2, T - t/2) e^{+i ω t}$

or

$x_W(τ, T) = x(T + t/2, T - t/2)$

of a 2-point function x. Returns a tuple of x_W and the corresponding axes (ω, T) or (τ, T), depending on the fourier keyword.

The motivation for the Wigner transformation is that, given an autocorrelation function x, it reduces to the spectral density function at all times T for stationary processes, yet it is fully equivalent to the non-stationary autocorrelation function. Therefore, the Wigner (distribution) function tells us, roughly, how the spectral density changes in time.

Optional keyword parameters

  • ts::AbstractVector: Time grid for x. Defaults to a UnitRange.
  • fourier::Bool: Whether to Fourier transform. Defaults to true.

Notes

The algorithm only works when ts – and consequently x – is equidistant.

References

https://en.wikipedia.org/wiki/Wigner_distribution_function

http://tftb.nongnu.org

source

Langreth's rules

KadanoffBaym.TimeOrderedGreenFunctionType
TimeOrderedGreenFunction(L::AbstractMatrix, G::AbstractMatrix)

A simple time-ordered Green function structure for a hassle-free computation of the Langreth rules.

Parameters

  • L::AbstractMatrix: The lesser component
  • G::AbstractMatrix: The greater component
source
KadanoffBaym.convFunction
conv(L::AbstractTimeOrderedGreenFunction, R::AbstractTimeOrderedGreenFunction, ws::UpperTriangular)

Calculates a time-convolution between time-ordered Green functions through the Langreth rules.

Parameters

  • L::AbstractTimeOrderedGreenFunction: The left time-ordered Green function
  • R::AbstractTimeOrderedGreenFunction: The right time-ordered Green function
  • ws::UpperTriangular: An upper-triangular weight matrix containing the integration weights
source

Citation

If you use KadanoffBaym.jl in your research, please cite our paper:

@Article{10.21468/SciPostPhysCore.5.2.030,
	title={{Adaptive Numerical Solution of Kadanoff-Baym Equations}},
	author={Francisco Meirinhos and Michael Kajan and Johann Kroha and Tim Bode},
	journal={SciPost Phys. Core},
	volume={5},
	issue={2},
	pages={30},
	year={2022},
	publisher={SciPost},
	doi={10.21468/SciPostPhysCore.5.2.030},
	url={https://scipost.org/10.21468/SciPostPhysCore.5.2.030},
}