SingleSiteDMRGEngine

Inheritance Diagram

Inheritance diagram of tenpy.algorithms.dmrg.SingleSiteDMRGEngine

Methods

SingleSiteDMRGEngine.__init__(psi, model, ...)

SingleSiteDMRGEngine.diag(theta_guess)

Diagonalize the effective Hamiltonian represented by self.

SingleSiteDMRGEngine.environment_sweeps(N_sweeps)

Perform N_sweeps sweeps without optimization to update the environment.

SingleSiteDMRGEngine.free_no_longer_needed_envs()

Remove no longer needed environments after an update.

SingleSiteDMRGEngine.get_resume_data([...])

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

SingleSiteDMRGEngine.get_sweep_schedule()

Define the schedule of the sweep.

SingleSiteDMRGEngine.init_env([model, ...])

(Re-)initialize the environment.

SingleSiteDMRGEngine.make_eff_H()

Create new instance of self.EffectiveH at self.i0 and set it to self.eff_H.

SingleSiteDMRGEngine.mixed_svd(theta)

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

SingleSiteDMRGEngine.mixer_activate()

Set self.mixer to the class specified by options['mixer'].

SingleSiteDMRGEngine.mixer_cleanup()

Cleanup the effects of a mixer.

SingleSiteDMRGEngine.plot_sweep_stats([...])

Plot sweep_stats to display the convergence with the sweeps.

SingleSiteDMRGEngine.plot_update_stats(axes)

Plot update_stats to display the convergence during the sweeps.

SingleSiteDMRGEngine.post_update_local(E0, ...)

Perform post-update actions.

SingleSiteDMRGEngine.prepare_svd(theta)

Transform theta into matrix for svd.

SingleSiteDMRGEngine.prepare_update()

Prepare self for calling update_local().

SingleSiteDMRGEngine.reset_stats([resume_data])

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

SingleSiteDMRGEngine.resume_run()

Resume a run that was interrupted.

SingleSiteDMRGEngine.run()

Run the DMRG simulation to find the ground state.

SingleSiteDMRGEngine.set_B(U, S, VH)

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

SingleSiteDMRGEngine.sweep([optimize, ...])

One 'sweep' of a the algorithm.

SingleSiteDMRGEngine.update_env(**update_data)

Update the left and right environments after an update of the state.

SingleSiteDMRGEngine.update_local(theta[, ...])

Perform site-update on the site i0.

Class Attributes and Properties

SingleSiteDMRGEngine.DMRG_params

SingleSiteDMRGEngine.engine_params

SingleSiteDMRGEngine.n_optimize

The number of sites to be optimized at once.

SingleSiteDMRGEngine.verbose

class tenpy.algorithms.dmrg.SingleSiteDMRGEngine(psi, model, options, **kwargs)[source]

Bases: tenpy.algorithms.dmrg.DMRGEngine

Engine for the single-site DMRG algorithm.

Parameters
  • psi (MPS) – Initial guess for the ground state, which is to be optimized in-place.

  • model (MPOModel) – The model representing the Hamiltonian for which we want to find the ground state.

  • options (dict) – Further optional parameters.

Options

config SingleSiteDMRGEngine
option summary

chi_list (from DMRGEngine) in DMRGEngine.reset_stats

A dictionary to gradually increase the `chi_max` parameter of [...]

chi_list_reactivates_mixer (from DMRGEngine) in DMRGEngine.sweep

If True, the mixer is reset/reactivated each time the bond dimension growth [...]

combine (from Sweep) in Sweep

Whether to combine legs into pipes. This combines the virtual and [...]

diag_method (from DMRGEngine) in DMRGEngine.run

Method to be used for diagonalzation, default ``'default'``. [...]

E_tol_max (from DMRGEngine) in DMRGEngine.run

See `E_tol_to_trunc`

E_tol_min (from DMRGEngine) in DMRGEngine.run

See `E_tol_to_trunc`

E_tol_to_trunc (from DMRGEngine) in DMRGEngine.run

It's reasonable to choose the Lanczos convergence criteria [...]

init_env_data (from Sweep) in DMRGEngine.init_env

Dictionary as returned by ``self.env.get_initialization_data()`` from [...]

lanczos_params (from Sweep) in Sweep

Lanczos parameters as described in :cfg:config:`Lanczos`.

max_E_err (from DMRGEngine) in DMRGEngine.run

Convergence if the change of the energy in each step [...]

max_hours (from DMRGEngine) in DMRGEngine.run

If the DMRG took longer (measured in wall-clock time), [...]

max_N_for_ED (from DMRGEngine) in DMRGEngine.diag

Maximum matrix dimension of the effective hamiltonian [...]

max_S_err (from DMRGEngine) in DMRGEngine.run

Convergence if the relative change of the entropy in each step [...]

max_sweeps (from DMRGEngine) in DMRGEngine.run

Maximum number of sweeps to be performed.

min_sweeps (from DMRGEngine) in DMRGEngine.run

Minimum number of sweeps to be performed. [...]

mixer (from DMRGEngine) in DMRGEngine.mixer_activate

Chooses the :class:`Mixer` to be used. [...]

mixer_params (from DMRGEngine) in DMRGEngine.mixer_activate

Mixer parameters as described in :cfg:config:`Mixer`.

N_sweeps_check (from DMRGEngine) in DMRGEngine.run

Number of sweeps to perform between checking convergence [...]

norm_tol (from DMRGEngine) in DMRGEngine.run

After the DMRG run, update the environment with at most [...]

norm_tol_iter (from DMRGEngine) in DMRGEngine.run

Perform at most `norm_tol_iter`*`update_env` sweeps to [...]

orthogonal_to (from Sweep) in DMRGEngine.init_env

Deprecated in favor of the `orthogonal_to` function argument (forwarded fro [...]

P_tol_max (from DMRGEngine) in DMRGEngine.run

See `P_tol_to_trunc`

P_tol_min (from DMRGEngine) in DMRGEngine.run

See `P_tol_to_trunc`

P_tol_to_trunc (from DMRGEngine) in DMRGEngine.run

It's reasonable to choose the Lanczos convergence criteria [...]

start_env (from Sweep) in DMRGEngine.init_env

Number of sweeps to be performed without optimization to update the environment.

sweep_0 (from DMRGEngine) in DMRGEngine.reset_stats

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

trunc_params (from Algorithm) in Algorithm

Truncation parameters as described in :cfg:config:`truncation`.

update_env (from DMRGEngine) in DMRGEngine.run

Number of sweeps without bond optimizaiton to update the [...]

chi_list

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.

Type

dict | None

eff_H

Effective two-site Hamiltonian.

Type

EffectiveH

mixer

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

Type

Mixer | None

shelve

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

Type

bool

sweeps

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

Type

int

time0

Time marker for the start of the run.

Type

float

update_stats

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

Type

dict

sweep_stats

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

Type

dict

EffectiveH[source]

alias of tenpy.algorithms.mps_common.OneSiteH

DefaultMixer[source]

alias of tenpy.algorithms.dmrg.SubspaceExpansion

prepare_svd(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(theta)[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 --
|                  |        |                   |         |

For a move to the left:

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

Note that theta lives on the same site i0 in both cases, but the sites of next_A and next_B depend on whether we move right or left. The returned U and VH have the same labels independent of that.

Without a mixer, this is done by a simple svd and truncation of Schmidt values of theta followed by the absorption of VH into next_B (U into next_A).

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

Parameters

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

Returns

  • U (Array) – Left-canonical part of theta. Labels '(vL.p)', 'vR'

  • S (1D 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.

  • VH (Array) – Right-canonical part of theta. Labels 'vL', '(p.vR)'.

  • err (TruncationError) – The truncation error introduced.

  • S_approx (ndarray) – Just the S if a 1D ndarray, or an approximation of the correct S (which was used for truncation) in case S is 2D Array.

set_B(U, S, VH)[source]

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

Parameters
  • U (Array) – Left and Right-canonical matrices as returned by the SVD.

  • VH (Array) – Left and Right-canonical matrices as returned by the SVD.

  • S (1D 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.

diag(theta_guess)[source]

Diagonalize the effective Hamiltonian represented by self.

option DMRGEngine.max_N_for_ED: int

Maximum matrix dimension of the effective hamiltonian up to which the 'default' diag_method uses ED instead of Lanczos.

option DMRGEngine.diag_method: str

One of the folloing strings:

‘default’

Same as 'lanczos' for large bond dimensions, but if the total dimension of the effective Hamiltonian does not exceed the DMRG parameter 'max_N_for_ED' it uses 'ED_block'.

‘lanczos’

lanczos() Default, the Lanczos implementation in 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_guess (Array) – Initial guess for the ground state of the effective Hamiltonian.

Returns

  • E0 (float) – Energy of the found ground state.

  • theta (Array) – Ground state of the effective Hamiltonian.

  • N (int) – Number of Lanczos iterations used. -1 if unknown.

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

environment_sweeps(N_sweeps)[source]

Perform N_sweeps sweeps without optimization to update the environment.

Parameters

N_sweeps (int) – Number of sweeps to run without optimization

free_no_longer_needed_envs()[source]

Remove no longer needed environments after an update.

This allows to minimize the number of environments to be kept. For large MPO bond dimensions, these environments are by far the biggest part in memory, so this is a valuable optimiztion to reduce memory requirements.

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 simulation from where we are now. It might contain an explicit copy of psi.

Return type

dict

get_sweep_schedule()[source]

Define the schedule of the sweep.

One ‘sweep’ is a full sequence from the leftmost site to the right and back.

Returns

schedule – 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 of the environments.

Return type

iterable of (int, bool, (bool, bool))

init_env(model=None, resume_data=None, orthogonal_to=None)[source]

(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
  • model (MPOModel) – The model representing the Hamiltonian for which we want to find the ground state. If None, keep the model used before.

  • resume_data (None | dict) – Given when resuming a simulation, as returned by get_resume_data(). Can contain another dict under the key init_env_data; the contents of init_env_data get passed as keyword arguments to the environment initializaiton.

  • orthogonal_to (None | list of MPS | list of dict) – List of other matrix product states to orthogonalize against. Instead of just the state, you can specify a dict with the state as ket and further keyword arguments for initializing the MPSEnvironment; the psi to be optimized is used as bra. Works only for finite or segment MPS; for infinite MPS it must be None. This can be used to find (a few) excited states as follows. First, run DMRG to find the ground state, and then run DMRG again while orthogonalizing against the ground state, which yields the first excited state (in the same symmetry sector), and so on. Note that resume_data['orthogonal_to'] takes precedence over the argument.

Options

Deprecated since version 0.6.0: Options LP, LP_age, RP and RP_age are now collected in a dictionary init_env_data with different keys init_LP, init_RP, age_LP, age_RP

Deprecated since version 0.8.0: Instead of passing the init_env_data as a option, it should be passed as dict entry of resume_data.

option Sweep.init_env_data: dict

Dictionary as returned by self.env.get_initialization_data() from get_initialization_data(). Deprecated, use the resume_data function/class argument instead.

option Sweep.orthogonal_to: list of MPS

Deprecated in favor of the orthogonal_to function argument (forwarded from the class argument) with the same effect.

option Sweep.start_env: int

Number of sweeps to be performed without optimization to update the environment.

Raises

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

make_eff_H()[source]

Create new instance of self.EffectiveH at self.i0 and set it to self.eff_H.

mixer_activate()[source]

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

option DMRGEngine.mixer: str | class | bool

Chooses the Mixer to be used. A string stands for one of the mixers defined in this module, a class is used as custom mixer. Default (None) uses no mixer, True uses DensityMatrixMixer for the 2-site case and SubspaceExpansion for the 1-site case. TwoSiteDMRGEngine only supports two-site mixers, but SingleSiteDMRGEngine supports both single-site and two-site mixers.

option DMRGEngine.mixer_params: dict

Mixer parameters as described in Mixer.

mixer_cleanup()[source]

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.

property n_optimize

The number of sites to be optimized at once.

Indirectly set by the class attribute EffectiveH and it’s length. For example, TwoSiteDMRGEngine uses the TwoSiteH and hence has n_optimize=2, while the SingleSiteDMRGEngine has n_optimize=1.

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

Plot sweep_stats to display the convergence with the sweeps.

Parameters
  • axes (matplotlib.axes.Axes) – The axes to plot into. Defaults to matplotlib.pyplot.gca()

  • xaxis (key of sweep_stats) – Key of sweep_stats to be used for the x-axis and y-axis of the plots.

  • yaxis (key of sweep_stats) – Key of sweep_stats to be used for the x-axis and y-axis of the plots.

  • y_exact (float) – 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(axes, xaxis='time', yaxis='E', y_exact=None, **kwargs)[source]

Plot update_stats to display the convergence during the sweeps.

Parameters
  • axes (matplotlib.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_exact (float) – 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(E0, age, N, ov_change, err, **update_data)[source]

Perform post-update actions.

Compute truncation energy and collect statistics.

Parameters

**update_data (dict) – What was returned by update_local().

prepare_update()[source]

Prepare self for calling update_local().

Returns

theta – Current best guess for the ground state, which is to be optimized. Labels are 'vL', 'p0', 'p1', 'vR', or combined versions of it (if self.combine). For single-site DMRG, the 'p1' label is missing.

Return type

Array

reset_stats(resume_data=None)[source]

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

option DMRGEngine.chi_list: dict | 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.

option DMRGEngine.sweep_0: int

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

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()[source]

Run the DMRG simulation to find the ground state.

Returns

  • E (float) – The energy of the resulting ground state MPS.

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

Options

option DMRGEngine.diag_method: str

Method to be used for diagonalzation, default 'default'. For possible arguments see DMRGEngine.diag().

option DMRGEngine.E_tol_to_trunc: float

It’s reasonable to choose the Lanczos convergence criteria 'E_tol' not many magnitudes lower than the current truncation error. Therefore, if E_tol_to_trunc is not None, we update E_tol of lanczos_params to max_E_trunc*E_tol_to_trunc, restricted to the interval [E_tol_min, E_tol_max], where max_E_trunc is the maximal energy difference due to truncation right after each Lanczos optimization during the sweeps.

option DMRGEngine.E_tol_max: float

See E_tol_to_trunc

option DMRGEngine.E_tol_min: float

See E_tol_to_trunc

option DMRGEngine.max_E_err: float

Convergence if the change of the energy in each step satisfies -Delta E / max(|E|, 1) < max_E_err. Note that this is also satisfied if Delta E > 0, i.e., if the energy increases (due to truncation).

option DMRGEngine.max_hours: float

If the DMRG took longer (measured in wall-clock time), ‘shelve’ the simulation, i.e. stop and return with the flag shelve=True.

option DMRGEngine.max_S_err: float

Convergence if the relative change of the entropy in each step satisfies |Delta S|/S < max_S_err

option DMRGEngine.max_sweeps: int

Maximum number of sweeps to be performed.

option DMRGEngine.min_sweeps: int

Minimum number of sweeps to be performed. Defaults to 1.5*N_sweeps_check.

option DMRGEngine.N_sweeps_check: int

Number of sweeps to perform between checking convergence criteria and giving a status update.

option DMRGEngine.norm_tol: float

After the DMRG run, update the environment with at most norm_tol_iter sweeps until np.linalg.norm(psi.norm_err()) < norm_tol.

option DMRGEngine.norm_tol_iter: float

Perform at most norm_tol_iter`*`update_env sweeps to converge the norm error below norm_tol. If the state is not converged after that, call canonical_form() instead.

option DMRGEngine.P_tol_to_trunc: float

It’s reasonable to choose the Lanczos convergence criteria 'P_tol' not many magnitudes lower than the current truncation error. Therefore, if P_tol_to_trunc is not None, we update P_tol of lanczos_params to max_trunc_err*P_tol_to_trunc, restricted to the interval [P_tol_min, P_tol_max], where max_trunc_err is the maximal truncation error (discarded weight of the Schmidt values) due to truncation right after each Lanczos optimization during the sweeps.

option DMRGEngine.P_tol_max: float

See P_tol_to_trunc

option DMRGEngine.P_tol_min: float

See P_tol_to_trunc

option DMRGEngine.update_env: int

Number of sweeps without bond optimizaiton to update the environment for infinite boundary conditions, performed every N_sweeps_check sweeps.

sweep(optimize=True, meas_E_trunc=False)[source]

One ‘sweep’ of a the algorithm.

Iteratate over the bond which is optimized, to the right and then back to the left to the starting point.

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

  • meas_E_trunc (bool, optional) – Whether to measure truncation energies.

Options

option DMRGEngine.chi_list_reactivates_mixer: bool

If True, the mixer is reset/reactivated each time the bond dimension growths due to DMRGEngine.chi_list.

Returns

  • max_trunc_err (float) – Maximal truncation error introduced.

  • max_E_trunc (None | float) – None if meas_E_trunc is False, else the maximal change of the energy due to the truncation.

update_env(**update_data)[source]

Update the left and right environments after an update of the state.

Parameters

**update_data – Whatever is returned by update_local().

update_local(theta, optimize=True)[source]

Perform site-update on the site i0.

Parameters
  • theta (Array) – Initial guess for the ground state of the effective Hamiltonian.

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

Returns

update_data – 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!

Return type

dict