PurificationTEBD

Inheritance Diagram

Inheritance diagram of tenpy.algorithms.purification.PurificationTEBD

Methods

PurificationTEBD.__init__(psi, model, ...)

PurificationTEBD.calc_U(order, delta_t[, ...])

see calc_U()

PurificationTEBD.disentangle(theta)

Disentangle theta before splitting with svd.

PurificationTEBD.disentangle_global([pair])

Try global disentangling by determining the maximally entangled pairs of sites.

PurificationTEBD.disentangle_global_nsite([n])

Perform a sweep through the system and disentangle with disentangle_n_site().

PurificationTEBD.disentangle_n_site(i, n, theta)

Generalization of disentangle() to n sites.

PurificationTEBD.estimate_RAM([...])

Gives an approximate prediction for the required memory usage.

PurificationTEBD.evolve(N_steps, dt)

Evolve by dt * N_steps.

PurificationTEBD.evolve_step(U_idx_dt, odd)

Updates either even or odd bonds in unit cell.

PurificationTEBD.get_resume_data([...])

Return necessary data to resume a run() interrupted at a checkpoint.

PurificationTEBD.prepare_evolve(dt)

Prepare an evolution step.

PurificationTEBD.resume_run()

Resume a run that was interrupted.

PurificationTEBD.run()

Perform a (real-)time evolution of psi by N_steps * dt.

PurificationTEBD.run_GS()

TEBD algorithm in imaginary time to find the ground state.

PurificationTEBD.run_evolution(N_steps, dt)

Perform a (real-)time evolution of psi by N_steps * dt.

PurificationTEBD.run_imaginary(beta)

Run imaginary time evolution to cool down by the given beta.

PurificationTEBD.suzuki_trotter_decomposition(...)

Returns list of necessary steps for the suzuki trotter decomposition.

PurificationTEBD.suzuki_trotter_time_steps(order)

Return time steps of U for the Suzuki Trotter decomposition of desired order.

PurificationTEBD.switch_engine(other_engine, *)

Initialize algorithm from another algorithm instance of a different class.

PurificationTEBD.update_bond(i, U_bond)

Updates the B matrices on a given bond.

PurificationTEBD.update_bond_imag(i, U_bond)

Update a bond with a (possibly non-unitary) U_bond.

PurificationTEBD.update_imag(N_steps[, ...])

Perform an update suitable for imaginary time evolution.

Class Attributes and Properties

PurificationTEBD.disent_iterations

For each bond the total number of iterations performed in any Disentangler.

PurificationTEBD.time_dependent_H

whether the algorithm supports time-dependent H

PurificationTEBD.trunc_err_bonds

truncation error introduced on each non-trivial bond.

class tenpy.algorithms.purification.PurificationTEBD(psi, model, options, **kwargs)[source]

Bases: TEBDEngine

Time evolving block decimation (TEBD) for purification MPS.

Parameters are the same as for Algorithm.

Options

config PurificationTEBD
option summary

delta_tau_list (from TEBDEngine) in PurificationTEBD.run_GS

A list of floats: the timesteps to be used. [...]

disentangle

Method for disentangling. [...]

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`.

option disentangle: None | str

Method for disentangling. See get_disentangler().

used_disentangler

The disentangler to be used on the auxiliary indices. Chosen by get_disentangler(), called with the TEBD parameter 'disentangle'. Defaults to the trivial disentangler for options['disentangle']=None.

Type:

Disentangler

_disent_iterations

Number of iterations performed on all bonds, including trivial bonds; length L.

Type:

1D ndarray

_guess_U_disent

Same index structure as self._U: for each two-site U of the physical time evolution the disentangler from the last application. Initialized to identities.

Type:

list of list of npc.Array

run_imaginary(beta)[source]

Run imaginary time evolution to cool down by the given beta.

Applies imaginary time evolution exp(-beta H) to psi.

Note that we don’t change the norm attribute of the MPS, i.e. normalization is preserved.

Parameters:

beta (float) – The inverse temperature beta = 1/T, by which we should cool down. We evolve to the closest multiple of options['dt'], see also evolved_time.

property disent_iterations

For each bond the total number of iterations performed in any Disentangler.

calc_U(order, delta_t, type_evo='real', E_offset=None)[source]

see calc_U()

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. This would look something like:

|           |             |
|     ... - B1  -  s  -  B2 - ...
|           |             |
|           |-------------|
|           |      U      |
|           |-------------|
|           |             |
Parameters:
  • i (int) – Bond index; we update the matrices at sites i-1, i.

  • U_bond (Array) – The bond operator which we apply to the wave function. We expect labels 'p0', 'p1', 'p0*', 'p1*' for U_bond.

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:
  • i (int) – Bond index; we update the matrices at sites i-1, i.

  • U_bond (Array) – The bond operator which we apply to the wave function. We expect labels 'p0', 'p1', 'p0*', 'p1*'.

Returns:

trunc_err – The error of the represented state which is introduced by the truncation during this update step.

Return type:

TruncationError

disentangle(theta)[source]

Disentangle theta before splitting with svd.

For the purification we write \(\rho_P = Tr_Q{|\psi_{P,Q}><\psi_{P,Q}|}\). Thus, we can actually apply any unitary to the auxiliary Q space of \(|\psi>\) without changing the result.

Note

We have to apply the same unitary to the ‘bra’ and ‘ket’ used for expectation values / correlation functions!

The behavior of this function is set by used_disentangler, which in turn is obtained from get_disentangler(options['disentangle']), see get_disentangler() for details on the syntax.

Parameters:

theta (Array) – Wave function to disentangle, with legs 'vL', 'vR', 'p0', 'p1', 'q0', 'q1'.

Returns:

  • theta_disentangled (Array) – Disentangled theta; npc.tensordot(U, theta, axes=[['q0*', 'q1*'], ['q0', 'q1']]).

  • U (Array) – The unitary used to disentangle theta, with labels 'q0', 'q1', 'q0*', 'q1*'. If no unitary was found/applied, it might also be None.

disentangle_global(pair=None)[source]

Try global disentangling by determining the maximally entangled pairs of sites.

Calculate the mutual information (in the auxiliary space) between two sites and determine where it is maximal. Disentangle these two sites with disentangle()

disentangle_global_nsite(n=2)[source]

Perform a sweep through the system and disentangle with disentangle_n_site().

Parameters:

n (int) – maximal number of sites to disentangle at once.

disentangle_n_site(i, n, theta)[source]

Generalization of disentangle() to n sites.

Simply group left and right n/2 physical legs, adjust labels, and apply disentangle() to disentangle the central bond. Recursively proceed to disentangle left and right parts afterwards. Scales (for even n) as \(O(\chi^3 d^n d^{n/2})\).

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:

float

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 with Us[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 save psi, model and options 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:

dict

prepare_evolve(dt)[source]

Prepare an evolution step.

This method is used to prepare repeated calls of evolve() given the model. For example, it may generate approximations of U=exp(-i H dt). To avoid overhead, it may cache the result depending on parameters/options; but it should always regenerate it if force_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 the run() 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 call run().

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 see TimeEvolutionAlgorithm.

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. Order 1 approximation is simply \(e^A a^B\). Order 2 is the “leapfrog” e^{A/2} e^B e^{A/2}. Order 4 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-steps d = suzuki_trotter_time_step(order) and the decomposition of H. They are chosen such that a subsequent application of exp(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:

list of (int, int)

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.

Parameters:

order (int) – The desired order of the Suzuki-Trotter decomposition.

Returns:

time_steps – We need U = exp(-i H_{even/odd} delta_t * dt) for the dt returned in this list.

Return type:

list of float

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 the TwoSiteDMRGEngine to the OneSiteDMRGEngine.

  • 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