quspin.tools.Floquet.Floquet

class quspin.tools.Floquet.Floquet(evo_dict, HF=False, UF=False, thetaF=False, VF=False, n_jobs=1, force_ONB=False)[source]

Bases: object

Calculates the Floquet spectrum, Floquet Hamiltonian and Floquet states.

Loops over the basis states to compute the Floquet unitary \(U_F\) (evolution operator over one period) for a periodically-driven system governed by the Hamiltonian \(H(t)=H(t+T)\):

\[U_F=U(T,0)=\mathcal{T}_t\exp\left(-i\int_0^T\mathrm{d}t H(t) \right)\]

with \(\mathcal{T}_t\exp\) denoting the time-ordered exponential.

Examples

Consider the following periodically driven spin-1/2 Hamiltonian

\[H(t) = \left\{ \begin{array}{cl} \sum_j J\sigma^z_{j+1}\sigma^z_j + h\sigma^z_j , & t\in[-T/4,T/4] \newline \sum_j g\sigma^x_j, & t \in[T/4,3T/4] \end{array} \right\} \mathrm{mod}\ T\]

where \(T=2\pi/\Omega\) is the drive period. We choose the starting point of the evolution (or equivalently – the driving phase) to be \(t=0\).

The following snippet of code shows how to calculate the Floquet eigenstates and the corresponding quasienergies, using evo_dict variable, case ii (see below).

 1from quspin.basis import spin_basis_1d  # Hilbert space spin basis
 2from quspin.tools.Floquet import Floquet, Floquet_t_vec  # Floquet Hamiltonian
 3import numpy as np  # generic math functions
 4
 5#
 6##### define model parameters #####
 7L = 10  # system size
 8J = 1.0  # spin interaction
 9g = 0.809  # transverse field
10h = 0.9045  # parallel field
11Omega = 4.5  # drive frequency
12
13
14#
15##### set up alternating Hamiltonians #####
16# define time-reversal symmetric periodic step drive
17def drive(t, Omega):
18    return np.sign(np.cos(Omega * t))
19
20
21drive_args = [Omega]
22# compute basis in the 0-total momentum and +1-parity sector
23basis = spin_basis_1d(L=L, a=1, kblock=0, pblock=1)
24# define PBC site-coupling lists for operators
25x_field_pos = [[+g, i] for i in range(L)]
26x_field_neg = [[-g, i] for i in range(L)]
27z_field = [[h, i] for i in range(L)]
28J_nn = [[J, i, (i + 1) % L] for i in range(L)]  # PBC
29# static and dynamic lists
30static = [["zz", J_nn], ["z", z_field], ["x", x_field_pos]]
31dynamic = [
32    ["zz", J_nn, drive, drive_args],
33    ["z", z_field, drive, drive_args],
34    ["x", x_field_neg, drive, drive_args],
35]
36# compute Hamiltonian
37H = 0.5 * hamiltonian(static, dynamic, dtype=np.float64, basis=basis)
38##### define time vector of stroboscopic times with 1 driving cycles and 10 points per cycle #####
39t = Floquet_t_vec(
40    Omega, 1, len_T=10
41)  # t.vals=times, t.i=initial time, t.T=drive period
42#
43##### calculate exact Floquet eigensystem #####
44t_list = (
45    np.array([0.0, t.T / 4.0, 3.0 * t.T / 4.0]) + np.finfo(float).eps
46)  # times to evaluate H
47dt_list = np.array(
48    [t.T / 4.0, t.T / 2.0, t.T / 4.0]
49)  # time step durations to apply H for
50Floq = Floquet(
51    {"H": H, "t_list": t_list, "dt_list": dt_list}, VF=True, force_ONB=True,
52)  # call Floquet class
53VF = Floq.VF  # read off Floquet states
54EF = Floq.EF  # read off quasienergies
__init__(evo_dict, HF=False, UF=False, thetaF=False, VF=False, n_jobs=1, force_ONB=False)[source]

Instantiates the Floquet class.

Parameters:
evo_dictdict

Dictionary which passes the different types of protocols to calculate the Floquet unitary. Depending on the protocol type, it contains the following keys:

  1. Periodic continuous protocol from a hamiltonian object.
    • H : hamiltonian object to generate the time evolution.

    • T : period of the protocol.

    • rtol : (optional) relative tolerance for the ODE solver. (default = 1E-9)

    • atol : (optional) absolute tolerance for the ODE solver. (default = 1E-9)

  2. Periodic step protocol from a hamiltonian object.
    • H : single hamiltonian object to generate the hamiltonians at each step. Periodic step drives can be encoded using a single function, e.g. \(\mathrm{sign}(\cos(\Omega t))\).

    • t_list : list of times to evaluate the hamiltonian at for each step.

    • dt_list : list of time step durations for each step of the evolution.

    • T: (optional) drive period used to compute the Floquet Hamiltonian H_F. If not specified, then T=sum(dt_list). Use this option for periodic delta kicks.

  3. Periodic step protocol from a list of hamiltonians.
    • H_list : list of matrices to evolve with.

    • dt_list : list of time step durations. Must be the same size as H_list.

    • T: (optional) drive period used to compute the Floquet Hamiltonian H_F. If not specified, then T=sum(dt_list). Use this option for periodic delta kicks.

HFbool

Set to True to calculate and return Floquet Hamiltonian under attribute _.HF. Default is False.

UFbool

Set to True to save evolution operator under attribute _.UF. Default is False.

thetaFbool

Set to True to save eigenvalues of the evolution operator (Floquet phases) under attribute _.thetaF. Default is False.

VFbool

Set to True to save Floquet states under attribute _.VF. Default is False.

n_jobsint, optional

Sets the number of processors which are used when looping over the basis states to compute the Floquet unitary. Default is False.

force_ONBbool

Set to True to run an extra QR decomposition to orthogonalize the Floquet states; only effective for VF=True.

Methods

__init__(evo_dict[, HF, UF, thetaF, VF, ...])

Instantiates the Floquet class.

Attributes

EF

ordered Floquet quasi-energies in interval \([-\Omega,\Omega]\).

HF

Floquet Hamiltonian.

T

drive period.

UF

Floquet unitary.

VF

Floquet eigenbasis (in columns).

thetaF

Floquet eigenphases.

property EF

ordered Floquet quasi-energies in interval \([-\Omega,\Omega]\).

Type:

numpy.ndarray(float)

property HF

Floquet Hamiltonian.

Requires __init__ argument HF=True.

Type:

numpy.ndarray(float)

property T

drive period.

Type:

float

property UF

Floquet unitary.

Requires __init__ argument UF=True.

Type:

numpy.ndarray(float)

property VF

Floquet eigenbasis (in columns).

Requires __init__ argument VF=True.

Type:

numpy.ndarray(float)

property thetaF

Floquet eigenphases.

Requires __init__ argument thetaF=True.

Type:

numpy.ndarray(float)