QRBasedTEBDEngine
full name: tenpy.algorithms.tebd.QRBasedTEBDEngine
parent module:
tenpy.algorithms.tebd
type: class
Inheritance Diagram
Methods
|
|
|
Calculate |
Gives an approximate prediction for the required memory usage. |
|
|
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.QRBasedTEBDEngine(psi, model, options, **kwargs)[source]
Bases:
TEBDEngine
Version of TEBD that relies on QR decompositions rather than SVD.
As introduced in arXiv:2212.09782.
Todo
To use use_eig_based_svd == True, which makes sense on GPU only, we need to implement the _eig_based_svd for “non-square” matrices. This means that \(M^{\dagger} M\) and \(M M^{\dagger}\) dont have the same size, and we need to disregard those eigenvectors of the larger one, that have eigenvalue zero, since we dont have corresponding eigenvalues of the smaller one.
Options
- config QRBasedTEBDEngine
option summary Expansion rate. The QR-based decomposition is carried out at an expanded bo [...]
Expansion rate at low ``chi``. [...]
Minimum bond dimension increase for each block. Default is `1`.
Whether the truncation error should be computed exactly. [...]
delta_tau_list (from TEBDEngine) 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.
E_offset (from TEBDEngine) in TEBDEngine
Energy offset to be applied in :meth:`calc_U`, see doc there. [...]
max_delta_t (from TEBDEngine) in TEBDEngine
Threshold for raising errors on too large time steps. Default ``1.0``. [...]
max_N_sites_per_ring (from Algorithm) in Algorithm
Threshold for raising errors on too many sites per ring. Default ``18``. [...]
max_trunc_err (from TimeEvolutionAlgorithm) in TimeDependentHAlgorithm.evolve
Threshold for raising errors on too large truncation errors. Default ``0.01 [...]
N_steps (from TEBDEngine) in PurificationTEBD.run_GS
Number of steps before measurement can be performed
order (from TEBDEngine) 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`.
start_trunc_err (from TimeEvolutionAlgorithm) in TimeEvolutionAlgorithm
Initial truncation error for :attr:`trunc_err`.
trunc_params (from Algorithm) in Algorithm
Truncation parameters as described in :cfg:config:`truncation`.
Whether the SVD of the bond matrix :math:`\Xi` should be carried out numeri [...]
- option cbe_expand: float
Expansion rate. The QR-based decomposition is carried out at an expanded bond dimension
eta = (1 + cbe_expand) * chi
, wherechi
is the bond dimension before the time step. Default is 0.1.
- option cbe_expand_0: float
Expansion rate at low
chi
. If given, the expansion rate decreases linearly fromcbe_expand_0
atchi == 1
tocbe_expand
atchi == trunc_params['chi_max']
, then remains constant. If not given, the expansion rate iscbe_expand
at allchi
.
- option cbe_min_block_increase: int
Minimum bond dimension increase for each block. Default is 1.
- option use_eig_based_svd: bool
Whether the SVD of the bond matrix \(\Xi\) should be carried out numerically via the eigensystem. This is faster on GPUs, but less accurate. It makes no sense to do this on CPU. It is currently not supported for update_imag. Default is False.
- option compute_err: bool
Whether the truncation error should be computed exactly. Compared to SVD-based TEBD, computing the truncation error is significantly more expensive. If True (default), the full error is computed. Otherwise, the truncation error is set to NaN.
- 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:
TruncationError
- 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:
TruncationError
- 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.
- estimate_RAM(mem_saving_factor=None)[source]
Gives an approximate prediction for the required memory usage.
This calculation is based on the requested bond dimension, the local Hilbert space dimension, the number of sites, and the boundary conditions.
- Parameters:
mem_saving_factor (float) – Represents the amount of RAM saved due to conservation laws. By default, it is ‘None’ and is extracted from the model automatically. However, this is only possible in a few cases and needs to be estimated in most cases. This is due to the fact that it is dependent on the model parameters. If one has a better estimate, one can pass the value directly. This value can be extracted by building the initial state psi (usually by performing DMRG) and then calling
print(psi.get_B(0).sparse_stats())
TeNPy will automatically print the fraction of nonzero entries in the first line, for example,6 of 16 entries (=0.375) nonzero
. This fraction corresponds to the mem_saving_factor; in our example, it is 0.375.- Returns:
usage – Required RAM in MB.
- Return type:
See also
tenpy.simulations.simulation.estimate_simulation_RAM
global function calling this.
- 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 evolution steps.
- Return type:
TruncationError
- 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:
TruncationError
- 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 algorithm run from where we are now. It might contain an explicit copy of 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.
- 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 behavior 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_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
- 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
.
- 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 optimized 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:
- 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.
- 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
- property trunc_err_bonds
truncation error introduced on each non-trivial bond.
- update_imag(N_steps, call_canonical_form=True)[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 able to preserve at least the orthonormality of the canoncial form.
- Parameters:
N_steps (int) – The number of steps for which the whole lattice should be updated.
call_canonical_from (bool) – The singular values saved in the MPS are not exactly correct after the update, since the non-unitary update on other bonds can change them. To fix this, we call psi.canonical_form at the end. Since this is about as a expensive as a single sweep, we allow to disable it, e.g. during the imaginary evolution looking for ground states where the intermediate results is not so critical.
- Returns:
trunc_err – The error of the represented state which is introduced due to the truncation during this sequence of update steps.
- Return type:
TruncationError