PurificationTEBD2

class tenpy.algorithms.purification_tebd.PurificationTEBD2(psi, model, TEBD_params)[source]

Bases: tenpy.algorithms.purification_tebd.PurificationTEBD

Similar as PurificationTEBD, but perform sweeps instead of brickwall.

Instead of the A-B pattern of even/odd bonds used in TEBD, perform sweeps similar as in DMRG for real-time evolution (similar as update_imag() does for imaginary time evolution).

Attributes
disent_iterations

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

trunc_err_bonds

truncation error introduced on each non-trivial bond.

Methods

calc_U(self, order, delta_t[, type_evo, …])

see calc_U()

disentangle(self, theta)

Disentangle theta before splitting with svd.

disentangle_global(self[, pair])

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

disentangle_global_nsite(self[, n])

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

disentangle_n_site(self, i, n, theta)

Generalization of disentangle() to n sites.

run(self)

(Real-)time evolution with TEBD (time evolving block decimation).

run_GS(self)

TEBD algorithm in imaginary time to find the ground state.

run_imaginary(self, beta)

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

suzuki_trotter_decomposition(order, N_steps)

Returns list of necessary steps for the suzuki trotter decomposition.

suzuki_trotter_time_steps(order)

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

update(self, N_steps)

Evolve by N_steps * U_param['dt'].

update_bond(self, i, U_bond)

Updates the B matrices on a given bond.

update_bond_imag(self, i, U_bond)

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

update_imag(self, N_steps)

Perform an update suitable for imaginary time evolution.

update_step(self, U_idx_dt, odd)

Updates bonds in unit cell.

update(self, N_steps)[source]

Evolve by N_steps * U_param['dt'].

Parameters
N_stepsint

The number of steps for which the whole lattice should be updated.

Returns
trunc_errTruncationError

The error of the represented state which is introduced due to the truncation during this sequence of update steps.

update_step(self, U_idx_dt, odd)[source]

Updates bonds in unit cell.

Depending on the choice of odd, perform a sweep to the left or right, updating once per site with a time step given by U_idx_dt.

Parameters
U_idx_dtint

Time step index in self._U, evolve with Us[i] = self.U[U_idx_dt][i] at bond (i-1,i).

oddbool/int

Indication of whether to update even (odd=False,0) or even (odd=True,1) sites

Returns
trunc_errTruncationError

The error of the represented state which is introduced due to the truncation during this sequence of update steps.

calc_U(self, order, delta_t, type_evo='real', E_offset=None)

see calc_U()

property disent_iterations

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

disentangle(self, theta)

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 auxiliar 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 behaviour of this function is set by used_disentangler, which in turn is obtained from get_disentangler(TEBD_params['disentangle']), see get_disentangler() for details on the syntax.

Parameters
thetaArray

Wave function to disentangle, with legs 'vL', 'vR', 'p0', 'p1', 'q0', 'q1'.

Returns
theta_disentangledArray

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

UArray

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(self, pair=None)

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

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

disentangle_global_nsite(self, n=2)

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(self, i, n, theta)

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})\).

run(self)

(Real-)time evolution with TEBD (time evolving block decimation).

The following (optional) parameters are read out from the TEBD_params.

key

type

description

dt

float

Time step.

order

int

Order of the algorithm.

The total error scales as O(t, dt^order).

N_steps

int

Number of time steps dt to evolve. (The Trotter decompositions of order > 1 are slightly more efficient if more than one step is performed at once.)

trunc_params

dict

Truncation parameters as described in truncate().

run_GS(self)

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.

The following (optional) parameters are read out from the TEBD_params. Use verbose=1 to print the used parameters during runtime.

key

type

description

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 the gradually decrease the time step.

order

int

Order of the Suzuki-Trotter decomposition.

N_steps

int

Number of steps before measurement can be performed

trunc_params

dict

Truncation parameters as described in truncate()

run_imaginary(self, beta)

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

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

Parameters
betafloat

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

static suzuki_trotter_decomposition(order, N_steps)

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
orderint

The desired order of the Suzuki-Trotter decomposition.

Returns
ST_decompositionlist of (int, int)

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

static suzuki_trotter_time_steps(order)

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

See suzuki_trotter_decomposition() for details.

Parameters
orderint

The desired order of the Suzuki-Trotter decomposition.

Returns
time_stepslist of float

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

property trunc_err_bonds

truncation error introduced on each non-trivial bond.

update_bond(self, i, U_bond)

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
iint

Bond index; we update the matrices at sites i-1, i.

U_bondArray

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

Returns
trunc_errTruncationError

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

update_bond_imag(self, i, U_bond)

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
iint

Bond index; we update the matrices at sites i-1, i.

U_bondArray

The bond operator which we apply to the wave function. We expect labels 'p0', 'p1', 'p0*', 'p1*'.

Returns
trunc_errTruncationError

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

update_imag(self, N_steps)

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_stepsint

The number of steps for which the whole lattice should be updated.

Returns
trunc_errTruncationError

The error of the represented state which is introduced due to the truncation during this sequence of update steps.