SingleSiteDMRGEngine

class tenpy.algorithms.dmrg.SingleSiteDMRGEngine(psi, model, engine_params)[source]

Bases: tenpy.algorithms.dmrg.DMRGEngine

‘Engine’ for the single-site DMRG algorithm.

Parameters
psiMPS

Initial guess for the ground state, which is to be optimized in-place.

modelMPOModel

The model representing the Hamiltonian for which we want to find the ground state.

engine_paramsdict

Further optional parameters. These are usually algorithm-specific, and thus should be described in subclasses.

Attributes
EffectiveHclass type

Class for the effective Hamiltonian (i.e., a subclass of EffectiveH. Has a length class attribute which specifies the number of sites updated at once (e.g., whether we do single-site vs. two-site DMRG).

chi_listdict | None

A dictionary to gradually increase the chi_max parameter of trunc_params. The key defines starting from which sweep chi_max is set to the value, e.g. {0: 50, 20: 100} uses chi_max=50 for the first 20 sweeps and chi_max=100 afterwards. Overwrites trunc_params[‘chi_list’]`. By default (None) this feature is disabled.

eff_HEffectiveH

Effective two-site Hamiltonian.

mixerMixer | None

If None, no mixer is used (anymore), otherwise the mixer instance.

shelvebool

If a simulation runs out of time (time.time() - start_time > max_seconds), the run will terminate with shelve = True.

sweepsint

The number of sweeps already performed. (Useful for re-start).

time0float

Time marker for the start of the run.

update_statsdict

A dictionary with detailed statistics of the convergence. For each key in the following table, the dictionary contains a list where one value is added each time Engine.update_bond() is called.

key

description

i0

An update was performed on sites i0, i0+1.

age

The number of physical sites involved in the simulation.

E_total

The total energy before truncation.

N_lanczos

Dimension of the Krylov space used in the lanczos diagonalization.

time

Wallclock time evolved since time0 (in seconds).

sweep_statsdict

A dictionary with detailed statistics of the convergence. For each key in the following table, the dictionary contains a list where one value is added each time Engine.sweep() is called (with optimize=True).

key

description

sweep

Number of sweeps performed so far.

E

The energy before truncation (as calculated by Lanczos).

S

Maximum entanglement entropy.

time

Wallclock time evolved since time0 (in seconds).

max_trunc_err

The maximum truncation error in the last sweep

max_E_trunc

Maximum change or Energy due to truncation in the last sweep.

max_chi

Maximum bond dimension used.

norm_err

Error of canonical form np.linalg.norm(psi.norm_test()).

Methods

diag(self, theta_guess)

Diagonalize the effective Hamiltonian represented by self.

environment_sweeps(self, N_sweeps)

Perform N_sweeps sweeps without optimization to update the environment.

get_sweep_schedule(self)

Define the schedule of the sweep.

init_env(self[, model])

(Re-)initialize the environment.

mixed_svd(self, theta, next_B)

Get (truncated) B from the new theta (as returned by diag).

mixer_activate(self)

Set self.mixer to the class specified by engine_params[‘mixer’].

mixer_cleanup(self)

Cleanup the effects of a mixer.

plot_sweep_stats(self[, axes, xaxis, yaxis, …])

Plot sweep_stats to display the convergence with the sweeps.

plot_update_stats(self, axes[, xaxis, …])

Plot update_stats to display the convergence during the sweeps.

post_update_local(self, update_data[, …])

Perform post-update actions.

prepare_svd(self, theta)

Transform theta into matrix for svd.

prepare_update(self)

Prepare self to represent the effective Hamiltonian on site i0.

reset_stats(self)

Reset the statistics, useful if you want to start a new sweep run.

run(self)

Run the DMRG simulation to find the ground state.

set_B(self, U, S, VH)

Update the MPS with the U, S, VH returned by self.mixed_svd.

sweep(self[, optimize, meas_E_trunc])

One ‘sweep’ of a sweeper algorithm.

update_LP(self, U)

Update left part of the environment.

update_RP(self, VH)

Update right part of the environment.

update_local(self, theta[, optimize, …])

Perform site-update on the site i0.

prepare_update(self)[source]

Prepare self to represent the effective Hamiltonian on site i0.

Returns
thetaArray

Current best guess for the ground state, which is to be optimized. Labels 'vL', 'p0', 'vR', or combined versions of it (if self.combine).

update_local(self, theta, optimize=True, meas_E_trunc=False)[source]

Perform site-update on the site i0.

Parameters
thetaArray

Initial guess for the ground state of the effective Hamiltonian.

optimizebool

Wheter we actually optimize to find the ground state of the effective Hamiltonian. (If False, just update the environments).

meas_E_truncbool

Wheter to measure the energy after truncation.

Returns
update_datadict

Data computed during the local update, as described in the following:

E0float

Total energy, obtained before truncation (if optimize=True), or after truncation (if optimize=False) (but never None).

Nint

Dimension of the Krylov space used for optimization in the lanczos algorithm. 0 if optimize=False.

ageint

Current size of the DMRG simulation: number of physical sites involved into the contraction.

U, VH: Array

U and VH returned by mixed_svd().

ov_change: float

Change in the wave function 1. - abs(<theta_guess|theta>) induced by diag(), not including the truncation!

prepare_svd(self, theta)[source]

Transform theta into matrix for svd.

In contrast with the 2-site engine, the matrix here depends on the direction we move, as we need ‘p’ to point away from the direction we are going in.

mixed_svd(self, theta, next_B)[source]

Get (truncated) B from the new theta (as returned by diag).

The goal is to split theta and truncate it. For a move to the right:

|   -- theta -- next_B --    ==>    -- U -- S -- VH -- next_B --
|        |      |                      |               |

For a move to the left:

|   -- next_B -- theta -- ==>    -- next_B -- U -- S -- VH --
|      |         |                  |                   |

The VH for right-move or U for left-move is absorebed into the next_B.

Without a mixer, this is done by a simple svd and truncation of Schmidt values of theta followed by the absorption of VH/U.

With a mixer, the state is perturbed before the SVD. The details of the perturbation are defined by the Mixer class.

Parameters
thetaArray

The optimized wave function, prepared for svd with prepare_svd(), i.e. with combined legs.

nextBArray

MPS tensor at the site that will be visited next.

Returns
UArray

Left-canonical part of theta. Labels '(vL.p0)', 'vR'.

S1D ndarray | 2D Array

Without mixer just the singluar values of the array; with mixer it might be a general matrix with labels 'vL', 'vR'; see comment above.

VHArray

Right-canonical part of theta. Labels 'vL', '(p0.vR)'.

errTruncationError

The truncation error introduced.

set_B(self, U, S, VH)[source]

Update the MPS with the U, S, VH returned by self.mixed_svd.

Parameters
U, VHArray

Left and Right-canonical matrices as returned by the SVD.

S1D array | 2D Array

The middle part returned by the SVD, theta = U S VH. Without a mixer just the singular values, with enabled mixer a 2D array.

mixer_activate(self)[source]

Set self.mixer to the class specified by engine_params[‘mixer’].

update_LP(self, U)[source]

Update left part of the environment.

The site at which to update the environment depends on the direction of the sweep. If we are sweeping right, update the invironment at i0+1. If we are sweeping left, update the environment at i0

Parameters
UArray

The U as returned by SVD, with combined legs, labels '(vL.p0)', 'vR' if self.move_right, else 'vL', '(p0.vR)'.

update_RP(self, VH)[source]

Update right part of the environment.

The site at which to update the environment depends on the direction of the sweep. If we are sweeping right, update the invironment at i0. If we are sweeping left, update the environment at i0-1

Parameters
VHArray

The VH as returned by SVD, with combined legs, labels '(vL.p0)', 'vR' if self.move_right, else 'vL', '(p0.vR)'.

diag(self, theta_guess)

Diagonalize the effective Hamiltonian represented by self.

The method used depends on the DMRG parameter diag_method.

diag_method

Function, comment

‘lanczos’

lanczos() Default, the Lanczos implementation of TeNPy

‘arpack’

lanczos_arpack() Based on scipy.linalg.sparse.eigsh(). Slower than ‘lanczos’, since it needs to convert the npc arrays to numpy arrays during each matvec, and possibly does many more iterations.

‘ED_block’

full_diag_effH() Contract the effective Hamiltonian to a (large!) matrix and diagonalize the block in the charge sector of the initial state. Preserves the charge sector of the explicitly conserved charges. However, if you don’t preserve a charge explicitly, it can break it. For example if you use a SpinChain({'conserve': 'parity'}), it could change the total “Sz”, but not the parity of ‘Sz’.

‘ED_all’

full_diag_effH() Contract the effective Hamiltonian to a (large!) matrix and diagonalize it completely. Allows to change the charge sector even for explicitly conserved charges. For example if you use a SpinChain({'conserve': 'Sz'}), it can change the total “Sz”.

Parameters
theta_guessArray

Initial guess for the ground state of the effective Hamiltonian.

Returns
E0float

Energy of the found ground state.

thetaArray

Ground state of the effective Hamiltonian.

Nint

Number of Lanczos iterations used. -1 if unknown.

ov_changefloat

Change in the wave function 1. - abs(<theta_guess|theta_diag>)

environment_sweeps(self, N_sweeps)

Perform N_sweeps sweeps without optimization to update the environment.

Parameters
N_sweepsint

Number of sweeps to run without optimization

get_sweep_schedule(self)

Define the schedule of the sweep.

One ‘sweep’ is a full sequence from the leftmost site to the right and back. Only those LP and RP that can be used later should be updated.

Returns
scheduleiterable of (int, bool, (bool, bool))

Schedule for the sweep. Each entry is (i0, move_right, (update_LP, update_RP)), where i0 is the leftmost of the self.EffectiveH.length sites to be updated in update_local(), move_right indicates whether the next i0 in the schedule is rigth (True) of the current one, and update_LP, update_RP indicate whether it is necessary to update the LP and RP. The latter are chosen such that the environment is growing for infinite systems, but we only keep the minimal number of environment tensors in memory.

init_env(self, model=None)

(Re-)initialize the environment.

This function is useful to (re-)start a Sweep with a slightly different model or different (engine) parameters. Note that we assume that we still have the same psi. Calls reset_stats().

Parameters
modelMPOModel

The model representing the Hamiltonian for which we want to find the ground state. If None, keep the model used before.

Raises
ValueError

If the engine is re-initialized with a new model, which legs are incompatible with those of hte old model.

mixer_cleanup(self)

Cleanup the effects of a mixer.

A sweep() with an enabled Mixer leaves the MPS psi with 2D arrays in S. To recover the originial form, this function simply performs one sweep with disabled mixer.

plot_sweep_stats(self, axes=None, xaxis='time', yaxis='E', y_exact=None, **kwargs)

Plot sweep_stats to display the convergence with the sweeps.

Parameters
axesmatplotlib.axes.Axes

The axes to plot into. Defaults to matplotlib.pyplot.gca()

xaxis, yaxiskey of sweep_stats

Key of sweep_stats to be used for the x-axis and y-axis of the plots.

y_exactfloat

Exact value for the quantity on the y-axis for comparison. If given, plot abs((y-y_exact)/y_exact) on a log-scale yaxis.

**kwargs :

Further keyword arguments given to axes.plot(...).

plot_update_stats(self, axes, xaxis='time', yaxis='E', y_exact=None, **kwargs)

Plot update_stats to display the convergence during the sweeps.

Parameters
axesmatplotlib.axes.Axes

The axes to plot into. Defaults to matplotlib.pyplot.gca()

xaxis'N_updates' | 'sweep' | keys of update_stats

Key of update_stats to be used for the x-axis of the plots. 'N_updates' is just enumerating the number of bond updates, and 'sweep' corresponds to the sweep number (including environment sweeps).

yaxis'E' | keys of update_stats

Key of update_stats to be used for the y-axisof the plots. For ‘E’, use the energy (per site for infinite systems).

y_exactfloat

Exact value for the quantity on the y-axis for comparison. If given, plot abs((y-y_exact)/y_exact) on a log-scale yaxis.

**kwargs :

Further keyword arguments given to axes.plot(...).

post_update_local(self, update_data, meas_E_trunc=False)

Perform post-update actions.

Compute truncation energy, remove LP/RP that are no longer needed and collect statistics.

Parameters
update_datadict

Data computed during the local update, as described in the following list.

meas_E_truncbool, optional

Wheter to measure the energy after truncation.

reset_stats(self)

Reset the statistics, useful if you want to start a new sweep run.

run(self)

Run the DMRG simulation to find the ground state.

Returns
Efloat

The energy of the resulting ground state MPS.

psiMPS

The MPS representing the ground state after the simluation, i.e. just a reference to psi.

sweep(self, optimize=True, meas_E_trunc=False)

One ‘sweep’ of a sweeper algorithm.

Iteratate over the bond which is optimized, to the right and then back to the left to the starting point. If optimize=False, don’t actually diagonalize the effective hamiltonian, but only update the environment.

Parameters
optimizebool, optional

Whether we actually optimize to find the ground state of the effective Hamiltonian. (If False, just update the environments).

meas_E_truncbool, optional

Whether to measure truncation energies.

Returns
max_trunc_errfloat

Maximal truncation error introduced.

max_E_truncNone | float

None if meas_E_trunc is False, else the maximal change of the energy due to the truncation.