DMRGEngine

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

Bases: tenpy.algorithms.mps_sweeps.Sweep

Generic ‘Engine’ for the single-site DMRG algorithm.

This engine is implemented as a subclass of Sweep. It contains all methods that are generic between SingleSiteDMRGEngine and TwoSiteDMRGEngine.

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 at local update-level. For each key in the following table, the dictionary contains a list where one value is added each time DMRGEngine.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).

ov_change

1. - abs(<theta_guess|theta_diag>), where |theta_guess> is the initial guess for the wave function and |theta_diag> is the untruncated wave function returned by diag().

sweep_statsdict

A dictionary with detailed statistics at the sweep level. 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 (excluding environment sweeps) performed so far.

N_updates

Number of updates (including environment 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.

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_update(self)

Prepare everything algorithm-specific to perform a local update.

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.

sweep(self[, optimize, meas_E_trunc])

One ‘sweep’ of a sweeper algorithm.

update_local(self, theta, \*\*kwargs)

Perform algorithm-specific local update.

run(self)[source]

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.

reset_stats(self)[source]

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

post_update_local(self, update_data, meas_E_trunc=False)[source]

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.

diag(self, theta_guess)[source]

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

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

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

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

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

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_activate(self)

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

It is expected that different algorithms have differen ways of implementing mixers (with different defaults). Thus, this is algorithm-specific.

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.

prepare_update(self)

Prepare everything algorithm-specific to perform a local update.

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.

update_local(self, theta, **kwargs)

Perform algorithm-specific local update.