quimb.tensor.tn1d.tebd¶
Tools for performing TEBD like algorithms in 1D.
Attributes¶
Classes¶
An simple interacting hamiltonian object used, for instance, in TEBD. |
|
Class implementing Time Evolving Block Decimation (TEBD) [1]. |
Functions¶
|
The out-of-time-ordered correlator (OTOC) generating by two local |
Module Contents¶
- class quimb.tensor.tn1d.tebd.LocalHam1D(L, H2, H1=None, cyclic=False)[source]¶
Bases:
quimb.tensor.tnag.tebd.LocalHamGenAn simple interacting hamiltonian object used, for instance, in TEBD. Once instantiated, the
LocalHam1Dhamiltonian stores a single term per pair of sites, cached versions of which can be retrieved likeH.get_gate_expm((i, i + 1), -1j * 0.5)etc.- Parameters:
L (int) – The size of the hamiltonian.
H2 (array_like or dict[tuple[int], array_like]) – The sum of interaction terms. If a dict is given, the keys should be nearest neighbours like
(10, 11), apart from any default term which should have the keyNone, and the values should be the sum of interaction terms for that interaction.H1 (array_like or dict[int, array_like], optional) – The sum of single site terms. If a dict is given, the keys should be integer sites, apart from any default term which should have the key
None, and the values should be the sum of single site terms for that site.cyclic (bool, optional) – Whether the hamiltonian has periodic boundary conditions or not.
- terms¶
The terms in the hamiltonian, combined from the inputs such that there is a single term per pair.
Examples
A simple, translationally invariant, interaction-only
LocalHam1D:>>> XX = pauli('X') & pauli('X') >>> YY = pauli('Y') & pauli('Y') >>> ham = LocalHam1D(L=100, H2=XX + YY)
The same, but with a translationally invariant field as well:
>>> Z = pauli('Z') >>> ham = LocalHam1D(L=100, H2=XX + YY, H1=Z)
Specifying a default interaction and field, with custom values set for some sites:
>>> H2 = {None: XX + YY, (49, 50): (XX + YY) / 2} >>> H1 = {None: Z, 49: 2 * Z, 50: 2 * Z} >>> ham = LocalHam1D(L=100, H2=H2, H1=H1)
Specifying the hamiltonian entirely through site specific interactions and fields:
>>> H2 = {(i, i + 1): XX + YY for i in range(99)} >>> H1 = {i: Z for i in range(100)} >>> ham = LocalHam1D(L=100, H2=H2, H1=H1)
See also
SpinHam1D- L¶
- cyclic = False¶
- build_mpo_propagator_trotterized(x, site_tag_id='I{}', tags=None, upper_ind_id='k{}', lower_ind_id='b{}', shape='lrud', contract_sites=True, **split_opts)[source]¶
Build an MPO representation of
expm(H * x), i.e. the imaginary or real time propagator of this local 1D hamiltonian, using a first order trotterized decomposition.- Parameters:
x (float) – The time to evolve for. Note this does not include the imaginary prefactor of the Schrodinger equation, so real
xcorresponds to imaginary time evolution, and vice versa.site_tag_id (str) – A string specifiying how to tag the tensors at each site. Should contain a
'{}'placeholder. It is used to generate the actual tags like:map(site_tag_id.format, range(len(arrays))).tags (str or sequence of str, optional) – Global tags to attach to all tensors.
upper_ind_id (str) – A string specifiying how to label the upper physical site indices. Should contain a
'{}'placeholder. It is used to generate the actual indices like:map(upper_ind_id.format, range(len(arrays))).lower_ind_id (str) – A string specifiying how to label the lower physical site indices. Should contain a
'{}'placeholder. It is used to generate the actual indices like:map(lower_ind_id.format, range(len(arrays))).shape (str, optional) – String specifying layout of the tensors. E.g. ‘lrud’ (the default) indicates the shape corresponds left-bond, right-bond, ‘up’ physical index, ‘down’ physical index. End tensors have either ‘l’ or ‘r’ dropped from the string if not periodic.
contract_sites (bool, optional) – Whether to contract all the decomposed factors at each site to yield a single tensor per site, by default True.
split_opts – Supplied to
tensor_split().
- class quimb.tensor.tn1d.tebd.TEBD(p0, H, dt=None, tol=None, t0=0.0, split_opts=None, progbar=True, imag=False)[source]¶
Class implementing Time Evolving Block Decimation (TEBD) [1].
[1] Guifré Vidal, Efficient Classical Simulation of Slightly Entangled Quantum Computations, PRL 91, 147902 (2003)
- Parameters:
p0 (MatrixProductState) – Initial state.
H (LocalHam1D or array_like) – Local terms or single dense hamiltonian representing the two body interaction. If so should have shape
(d * d, d * d), wheredis the physical dimension ofp0.dt (float, optional) – Default time step, cannot be set as well as
tol.tol (float, optional) – Default target error for each evolution, cannot be set as well as
dt, which will instead be calculated from the trotter order, length of time, and hamiltonian norm.t0 (float, optional) – Initial time. Defaults to 0.0.
split_opts (dict, optional) – Compression options applied for splitting after gate application, see
tensor_split().imag (bool, optional) – Enable imaginary time evolution. Defaults to
False.
See also
- _pt¶
- L¶
- H¶
- cyclic¶
- _ham_norm¶
- _err = 0.0¶
- tol = None¶
- imag = False¶
- progbar = True¶
- split_opts¶
- property pt¶
The MPS state of the system at the current time.
- property err¶
- choose_time_step(tol, T, order)[source]¶
Trotter error is
~ (T / dt) * dt^(order + 1). Invert to find desired time step, and scale by norm of interaction term.
- _get_gate_from_ham(dt_frac, sites)[source]¶
Get the unitary (exponentiated) gate for fraction of timestep
dt_fracand sitessites, cached.
- sweep(direction, dt_frac, dt=None, queue=False)[source]¶
Perform a single sweep of gates and compression. This shifts the orthonognality centre along with the gates as they are applied and split.
- TARGET_TOL = 1e-13¶
- update_to(T, dt=None, tol=None, order=4, progbar=None)[source]¶
Update the state to time
T.- Parameters:
- at_times(ts, dt=None, tol=None, order=4, progbar=None)[source]¶
Generate the time evolved state at each time in
ts.- Parameters:
ts (sequence of float) – The times to evolve to and yield the state at.
dt (float, optional) – Time step to use. Can’t be set as well as
tol.tol (float, optional) – Tolerance for whole evolution. Can’t be set as well as
dt.order (int, optional) – Trotter order to use.
progbar (bool, optional) – Manually turn the progress bar off.
- Yields:
pt (MatrixProductState) – The state at each of the times in
ts. This is a copy of internal state used, so inplace changes can be made to it.
- quimb.tensor.tn1d.tebd.OTOC_local(psi0, H, H_back, ts, i, A, j=None, B=None, initial_eigenstate='check', **tebd_opts)[source]¶
The out-of-time-ordered correlator (OTOC) generating by two local operator A and B acting on site ‘i’, note it’s a function of time.
- Parameters:
psi0 (MatrixProductState) – The initial state in MPS form.
H (LocalHam1D) – The Hamiltonian for forward time-evolution.
H_back (LocalHam1D) – The Hamiltonian for backward time-evolution, should have only sign difference with ‘H’.
ts (sequence of float) – The time to evolve to.
i (int) – The site where the local operators acting on.
A (array) – The operator to act with.
initial_eigenstate ({'check', Flase, True}) – To check the psi0 is or not eigenstate of operator B. If psi0 is the eigenstate of B, it will run a simpler version of OTOC calculation automatically.
- Return type: