RandomUnitaryEvolution

  • full name: tenpy.algorithms.tebd.RandomUnitaryEvolution

  • parent module: tenpy.algorithms.tebd

  • type: class

Inheritance Diagram

Inheritance diagram of tenpy.algorithms.tebd.RandomUnitaryEvolution

Methods

RandomUnitaryEvolution.__init__(psi, options)

Initialize self.

RandomUnitaryEvolution.calc_U()

Draw new random two-site unitaries replacing the usual U of TEBD.

RandomUnitaryEvolution.get_resume_data()

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

RandomUnitaryEvolution.resume_run()

Resume a run that was interrupted.

RandomUnitaryEvolution.run()

Time evolution with TEBD and random two-site unitaries.

RandomUnitaryEvolution.run_GS()

TEBD algorithm in imaginary time to find the ground state.

RandomUnitaryEvolution.suzuki_trotter_decomposition(…)

Returns list of necessary steps for the suzuki trotter decomposition.

RandomUnitaryEvolution.suzuki_trotter_time_steps(order)

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

RandomUnitaryEvolution.update(N_steps)

Apply N_steps random two-site unitaries to each bond (in even-odd pattern).

RandomUnitaryEvolution.update_bond(i, U_bond)

Updates the B matrices on a given bond.

RandomUnitaryEvolution.update_bond_imag(i, …)

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

RandomUnitaryEvolution.update_imag(N_steps)

Perform an update suitable for imaginary time evolution.

RandomUnitaryEvolution.update_step(U_idx_dt, odd)

Updates either even or odd bonds in unit cell.

Class Attributes and Properties

RandomUnitaryEvolution.TEBD_params

RandomUnitaryEvolution.trunc_err_bonds

truncation error introduced on each non-trivial bond.

RandomUnitaryEvolution.verbose

class tenpy.algorithms.tebd.RandomUnitaryEvolution(psi, options)[source]

Bases: tenpy.algorithms.tebd.TEBDEngine

Evolution of an MPS with random two-site unitaries in a TEBD-like fashion.

Instead of using a model Hamiltonian, this TEBD engine evolves with random two-site unitaries. These unitaries are drawn according to the Haar measure on unitaries obeying the conservation laws dictated by the conserved charges. If no charge is preserved, this distribution is called circular unitary ensemble (CUE), see CUE().

On one hand, such an evolution is of interest in recent research (see eg. arXiv:1710.09827). On the other hand, it also comes in handy to “randomize” an initial state, e.g. for DMRG. Note that the entanglement grows very quickly, choose the truncation paramters accordingly!

Parameters
  • psi (MPS) – Initial state to be time evolved. Modified in place.

  • options (dict) – See below for details.

Options

config RandomUnitaryEvolution
option summary

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.

N_steps

Number of two-site unitaries to be applied on each bond.

order (from TEBDEngine) in PurificationTEBD.run_GS

Order of the Suzuki-Trotter decomposition.

start_time (from TimeEvolutionAlgorithm) in TimeEvolutionAlgorithm

Initial value for :attr:`evolved_time`.

start_trunc_err (from TEBDEngine) in TEBDEngine

Initial truncation error for :attr:`trunc_err`.

trunc_params

Truncation parameters as described in :cfg:config:`truncate`

option N_steps: int

Number of two-site unitaries to be applied on each bond.

option trunc_params: dict

Truncation parameters as described in truncate

Examples

One can initialize a “random” state with total Sz = L//2 as follows:

>>> from tenpy.algorithms.tebd import RandomUnitaryEvolution
>>> from tenpy.networks.mps import MPS
>>> L = 8
>>> spin_half = tenpy.networks.site.SpinHalfSite(conserve='Sz')
>>> psi = MPS.from_product_state([spin_half]*L, ["up", "down"]*(L//2), bc='finite')  # Neel
>>> print(psi.chi)
[1, 1, 1, 1, 1, 1, 1]
>>> options = dict(N_steps=2, trunc_params={'chi_max':10})
>>> eng = RandomUnitaryEvolution(psi, options)
>>> eng.run()
>>> print(psi.chi)
[2, 4, 8, 10, 8, 4, 2]
>>> psi.canonical_form()  # a good idea if there was a truncation necessary.

The “random” unitaries preserve the specified charges, e.g. here we have Sz-conservation. If you start in a sector of all up spins, the random unitaries can only apply a phase:

>>> psi2 = MPS.from_product_state([spin_half]*L, [0]*L, bc='finite')  # all spins up
>>> print(psi2.chi)
[1, 1, 1, 1, 1, 1, 1]
>>> eng2 = RandomUnitaryEvolution(psi2, options)
>>> eng2.run()  # random unitaries respect Sz conservation -> we stay in all-up sector
>>> print(psi2.chi)  # still a product state, not really random!!!
[1, 1, 1, 1, 1, 1, 1]
run()[source]

Time evolution with TEBD and random two-site unitaries.

calc_U()[source]

Draw new random two-site unitaries replacing the usual U of TEBD.

update(N_steps)[source]

Apply N_steps random two-site unitaries to each bond (in even-odd pattern).

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

TruncationError

get_resume_data()[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(resume_data)

An algorithm which doesn’t support this should override resume_run to raise an Error.

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.

Return type

dict

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 behaviour implemented here is to just call run().

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 the 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_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 optmized 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

property trunc_err_bonds

truncation error introduced on each non-trivial bond.

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 correponding tensor networks look like this:

|           --S--B1--B2--           --B1--B2--
|                |   |                |   |
|     theta:     U_bond        C:     U_bond
|                |   |                |   |
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

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

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

TruncationError

update_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