Time evolution
Schwinger.jl supports time-evolving states using all three backends:
- ED: Uses Krylov methods from
KrylovKit.jl - ITensors: Uses the time-dependent variational principle (TDVP) from
ITensorMPS.jl - MPSKit: Uses TDVP from
MPSKit.jl
The evolve function evolves a state forwards in time, and monitors any given observables. It returns a final state along with a DataFrame of the observables. For example, here is a simulation of flux unwinding.
using Schwinger, Plots
lat = Lattice(10; F = 1, L = 2, periodic = true)
gs = groundstate(Hamiltonian(lat; backend=:ED))
_, df = evolve(WilsonLoop(lat; backend=:ED) * gs, 15;
nsteps = 30,
observable = (ψ, t) -> sum(electricfields(ψ))/10
)
scatter(df.time, df.observable, xlabel = "gt", ylabel = "Average electric field", label = "Schwinger.jl")
plot!(0:.1:15, [cos(t/√(π)) for t in 0:.1:15], label = "Exact")Schwinger.evolve — Function
evolve(state::EDState, t::Real; nsteps::Int = 1, tol = 1E-12, observable = nothing, kwargs...)Evolve an exact diagonalization state by time t using matrix exponentiation.
Arguments
state::EDState: The state to evolve.t::Real: The time to evolve by.nsteps::Int = 1: Number of time steps to divide the evolution into.tol = 1E-12: Tolerance for the matrix exponentiation.observable::Union{Nothing,Function,Dict} = nothing: Observable(s) to monitor during evolution. Can be a single function(state, t) -> valueor a dictionary of name => function pairs.
Returns
EDState: The evolved state.obs: Observer object containing the history of observables and metadata.
Examples
evolved_state, obs = evolve(state, 1.0; nsteps=10, observable=s->energy(s))evolve(state::ITensorState, t::Real; nsteps::Int = 1, observable = nothing, kwargs...)Evolve an ITensor MPS state by imaginary time t using TDVP.
Arguments
state::ITensorState: The state to evolve.t::Real: The time to evolve by (will be multiplied by -im for imaginary time evolution).nsteps::Int = 1: Number of time steps for the TDVP algorithm.observable::Union{Nothing,Function,Dict} = nothing: Observable(s) to monitor during evolution. Can be a single function(state, t) -> valueor a dictionary of name => function pairs.kwargs...: Additional keyword arguments passed to the TDVP algorithm.
Returns
ITensorState: The evolved state.obs: Observer object containing the history of observables and metadata.
Notes
- Uses ITensors.jl TDVP algorithm.
- Time is tracked as real values (divided by -im) in the observer.
evolve(state::MPSKitState, t::Real; nsteps::Int = 1, dt = nothing, observable = nothing, kwargs...)Evolve an MPSKit MPS state by time t using MPSKit's TDVP algorithm.
Arguments
state::MPSKitState: The state to evolve.t::Real: The time to evolve by.nsteps::Int = 1: Number of time steps to divide the evolution into.dt::Union{Nothing,Real} = nothing: Time step size (currently unused, evolution uses t/nsteps).observable::Union{Nothing,Function,Dict} = nothing: Observable(s) to monitor during evolution. Can be a single function(state, t) -> valueor a dictionary of name => function pairs.kwargs...: Additional keyword arguments.
Returns
MPSKitState: The evolved state.obs: Observer object containing the history of observables and metadata.
Notes
- Uses MPSKit.jl TDVP algorithm.
- Evolves with imaginary time (-im*t/nsteps).