TEBDEngine¶
full name: tenpy.algorithms.tebd.TEBDEngine
parent module:
tenpy.algorithms.tebd
type: class
Inheritance Diagram

Methods
|
|
|
Calculate |
|
Evolve by |
|
Updates either even or odd bonds in unit cell. |
|
Return necessary data to resume a |
Prepare an evolution step. |
|
Resume a run that was interrupted. |
|
Perform a (real-)time evolution of |
|
TEBD algorithm in imaginary time to find the ground state. |
|
|
Perform a (real-)time evolution of |
Returns list of necessary steps for the suzuki trotter decomposition. |
|
Return time steps of U for the Suzuki Trotter decomposition of desired order. |
|
|
Initialize algorithm from another algorithm instance of a different class. |
|
Updates the B matrices on a given bond. |
|
Update a bond with a (possibly non-unitary) U_bond. |
|
Perform an update suitable for imaginary time evolution. |
Class Attributes and Properties
|
|
whether the algorithm supports time-dependent H |
|
truncation error introduced on each non-trivial bond. |
|
|
- class tenpy.algorithms.tebd.TEBDEngine(psi, model, options, **kwargs)[source]¶
Bases:
TimeEvolutionAlgorithm
Time Evolving Block Decimation (TEBD) algorithm.
Parameters are the same as for
Algorithm
.Deprecated since version 0.6.0: Renamed parameter/attribute TEBD_params to
options
.Options
- config TEBDEngine¶
option summary delta_tau_list in PurificationTEBD.run_GS
A list of floats: the timesteps to be used. [...]
dt (from TimeEvolutionAlgorithm) in TimeEvolutionAlgorithm
Minimal time step by which to evolve.
Energy offset to be applied in :meth:`calc_U`, see doc there. [...]
N_steps in PurificationTEBD.run_GS
Number of steps before measurement can be performed
order in PurificationTEBD.run_GS
Order of the Suzuki-Trotter decomposition.
preserve_norm (from TimeEvolutionAlgorithm) in TimeEvolutionAlgorithm
Whether the state will be normalized to its initial norm after each time st [...]
start_time (from TimeEvolutionAlgorithm) in TimeEvolutionAlgorithm
Initial value for :attr:`evolved_time`.
Initial truncation error for :attr:`trunc_err`.
trunc_params (from Algorithm) in Algorithm
Truncation parameters as described in :cfg:config:`truncation`.
-
option start_trunc_err:
TruncationError
¶ Initial truncation error for
trunc_err
.
- option order: int¶
Order of the algorithm. The total error for evolution up to a fixed time t scales as
O(t*dt^order)
.
-
option start_trunc_err:
- trunc_err¶
The error of the represented state which is introduced due to the truncation during the sequence of update steps.
- Type
- model¶
The model defining the Hamiltonian.
- Type
- _U¶
Exponentiated H_bond (bond Hamiltonians), i.e. roughly
exp(-i H_bond dt_i)
. First list for different dt_i as necessary for the chosen order, second list for the L different bonds.- Type
list of list of
Array
- _U_param¶
A dictionary containing the information of the latest created _U. We don’t recalculate _U if those parameters didn’t change.
- Type
- _trunc_err_bonds¶
The local truncation error introduced at each bond, ignoring the errors at other bonds. The i-th entry is left of site i.
- Type
list of
TruncationError
- _update_index¶
The indices
i_dt,i_bond
ofU_bond = self._U[i_dt][i_bond]
during update_step.
- property trunc_err_bonds¶
truncation error introduced on each non-trivial bond.
- run_GS()[source]¶
TEBD algorithm in imaginary time to find the ground state.
Note
It is almost always more efficient (and hence advisable) to use DMRG. This algorithms can nonetheless be used quite well as a benchmark and for comparison.
- option TEBDEngine.delta_tau_list: list¶
A list of floats: the timesteps to be used. Choosing a large timestep delta_tau introduces large (Trotter) errors, but a too small time step requires a lot of steps to reach
exp(-tau H) --> |psi0><psi0|
. Therefore, we start with fairly large time steps for a quick time evolution until convergence, and then gradually decrease the time step.
- option TEBDEngine.order: int¶
Order of the Suzuki-Trotter decomposition.
- option TEBDEngine.N_steps: int¶
Number of steps before measurement can be performed
- static suzuki_trotter_time_steps(order)[source]¶
Return time steps of U for the Suzuki Trotter decomposition of desired order.
See
suzuki_trotter_decomposition()
for details.
- static suzuki_trotter_decomposition(order, N_steps)[source]¶
Returns list of necessary steps for the suzuki trotter decomposition.
We split the Hamiltonian as \(H = H_{even} + H_{odd} = H[0] + H[1]\). The Suzuki-Trotter decomposition is an approximation \(\exp(t H) \approx prod_{(j, k) \in ST} \exp(d[j] t H[k]) + O(t^{order+1 })\).
- Parameters
order (
1, 2, 4, '4_opt'
) – The desired order of the Suzuki-Trotter decomposition. Order1
approximation is simply \(e^A a^B\). Order2
is the “leapfrog” e^{A/2} e^B e^{A/2}. Order4
is the fourth-order from [suzuki1991] (also referenced in [schollwoeck2011]), and'4_opt'
gives the optmized version of Equ. (30a) in [barthel2020].- Returns
ST_decomposition – Indices
j, k
of the time-stepsd = suzuki_trotter_time_step(order)
and the decomposition of H. They are chosen such that a subsequent application ofexp(d[j] t H[k])
to a given state|psi>
yields(exp(N_steps t H[k]) + O(N_steps t^{order+1}))|psi>
.- Return type
- prepare_evolve(dt)[source]¶
Prepare an evolution step.
This method is used to prepare repeated calls of
evolve()
given themodel
. For example, it may generate approximations ofU=exp(-i H dt)
. To avoid overhead, it may cache the result depending on parameters/options; but it should always regenerate it ifforce_prepare_evolve
is set.- Parameters
dt (float) – The time step to be used.
- calc_U(order, delta_t, type_evo='real', E_offset=None)[source]¶
Calculate
self.U_bond
fromself.model.H_bond
.This function calculates
U_bond = exp(-i dt (H_bond-E_offset_bond))
fortype_evo='real'
, orU_bond = exp(- dt H_bond)
fortype_evo='imag'
.
For first order (in delta_t), we need just one
dt=delta_t
. Higher order requires smaller dt steps, as given bysuzuki_trotter_time_steps()
.- Parameters
order (int) – Trotter order calculated U_bond. See update for more information.
delta_t (float) – Size of the time-step used in calculating U_bond
type_evo (
'imag' | 'real'
) – Determines whether we perform real or imaginary time-evolution.E_offset (None | list of float) – Possible offset added to H_bond for real-time evolution.
- evolve(N_steps, dt)[source]¶
Evolve by
dt * N_steps
.- Parameters
N_steps (int) – The number of steps for which the whole lattice should be updated.
dt (float) – The time step; but really this was already used in
prepare_evolve()
.
- Returns
trunc_err – The error of the represented state which is introduced due to the truncation during this sequence of evolvution steps.
- Return type
- evolve_step(U_idx_dt, odd)[source]¶
Updates either even or odd bonds in unit cell.
Depending on the choice of p, this function updates all even (
E
, odd=False,0) or odd (O
) (odd=True,1) bonds:| - B0 - B1 - B2 - B3 - B4 - B5 - B6 - | | | | | | | | | | |----| |----| |----| | | | E | | E | | E | | | |----| |----| |----| | |----| |----| |----| | | | O | | O | | O | | | |----| |----| |----| |
Note that finite boundary conditions are taken care of by having
Us[0] = None
.- Parameters
U_idx_dt (int) – Time step index in
self._U
, evolve withUs[i] = self.U[U_idx_dt][i]
at bond(i-1,i)
.odd (bool/int) – Indication of whether to update even (
odd=False,0
) or even (odd=True,1
) sites
- Returns
trunc_err – The error of the represented state which is introduced due to the truncation during this sequence of update steps.
- Return type
- update_bond(i, U_bond)[source]¶
Updates the B matrices on a given bond.
Function that updates the B matrices, the bond matrix s between and the bond dimension chi for bond i. The corresponding tensor networks look like this:
| --S--B1--B2-- --B1--B2-- | | | | | | theta: U_bond C: U_bond | | | | |
- Parameters
- Returns
trunc_err – The error of the represented state which is introduced by the truncation during this update step.
- Return type
- update_imag(N_steps)[source]¶
Perform an update suitable for imaginary time evolution.
Instead of the even/odd brick structure used for ordinary TEBD, we ‘sweep’ from left to right and right to left, similar as DMRG. Thanks to that, we are actually able to preserve the canonical form.
- Parameters
N_steps (int) – The number of steps for which the whole lattice should be updated.
- Returns
trunc_err – The error of the represented state which is introduced due to the truncation during this sequence of update steps.
- Return type
- update_bond_imag(i, U_bond)[source]¶
Update a bond with a (possibly non-unitary) U_bond.
Similar as
update_bond()
; but after the SVD just keep the A, S, B canonical form. In that way, one can sweep left or right without using old singular values, thus preserving the canonical form during imaginary time evolution.- Parameters
- Returns
trunc_err – The error of the represented state which is introduced by the truncation during this update step.
- Return type
- get_resume_data(sequential_simulations=False)[source]¶
Return necessary data to resume a
run()
interrupted at a checkpoint.At a
checkpoint
, you can savepsi
,model
andoptions
along with the data returned by this function. When the simulation aborts, you can resume it using this saved data with:eng = AlgorithmClass(psi, model, options, resume_data=resume_data) eng.resume_run()
An algorithm which doesn’t support this should override resume_run to raise an Error.
- Parameters
sequential_simulations (bool) – If True, return only the data for re-initializing a sequential simulation run, where we “adiabatically” follow the evolution of a ground state (for variational algorithms), or do series of quenches (for time evolution algorithms); see
run_seq_simulations()
.- Returns
resume_data – Dictionary with necessary data (apart from copies of psi, model, options) that allows to continue the simulation from where we are now. It might contain an explicit copy of psi.
- Return type
- resume_run()[source]¶
Resume a run that was interrupted.
In case we saved an intermediate result at a
checkpoint
, this function allows to resume therun()
of the algorithm (after re-initialization with the resume_data). Since most algorithms just have a while loop with break conditions, the default behaviour implemented here is to just callrun()
.
- run()[source]¶
Perform a (real-)time evolution of
psi
by N_steps * dt.You probably want to call this in a loop along with measurements. The recommended way to do this is via the
RealTimeEvolution
.
- run_evolution(N_steps, dt)[source]¶
Perform a (real-)time evolution of
psi
by N_steps * dt.This is the inner part of
run()
without the logging. For parameters seeTimeEvolutionAlgorithm
.
- classmethod switch_engine(other_engine, *, options=None, **kwargs)[source]¶
Initialize algorithm from another algorithm instance of a different class.
You can initialize one engine from another, not too different subclasses. Internally, this function calls
get_resume_data()
to extract data from the other_engine and then initializes the new class.Note that it transfers the data without making copies in most case; even the options! Thus, when you call run() on one of the two algorithm instances, it will modify the state, environment, etc. in the other. We recommend to make the switch as
engine = OtherSubClass.switch_engine(engine)
directly replacing the reference.- Parameters
cls (class) – Subclass of
Algorithm
to be initialized.other_engine (
Algorithm
) – The engine from which data should be transferred. Another, but not too different algorithm subclass-class; e.g. you can switch from theTwoSiteDMRGEngine
to theOneSiteDMRGEngine
.options (None | dict-like) – If not None, these options are used for the new initialization. If None, take the options from the other_engine.
**kwargs – Further keyword arguments for class initialization. If not defined, resume_data is collected with
get_resume_data()
.
- time_dependent_H = False¶
whether the algorithm supports time-dependent H