EngineCombine

Inheritance Diagram

Inheritance diagram of tenpy.algorithms.dmrg.EngineCombine

Methods

EngineCombine.__init__(psi, model, DMRG_params)

Initialize self.

EngineCombine.diag(theta_guess)

Diagonalize the effective Hamiltonian represented by self.

EngineCombine.environment_sweeps(N_sweeps)

Perform N_sweeps sweeps without optimization to update the environment.

EngineCombine.get_sweep_schedule()

Define the schedule of the sweep.

EngineCombine.init_env([model])

(Re-)initialize the environment.

EngineCombine.make_eff_H()

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

EngineCombine.mixed_svd(theta)

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

EngineCombine.mixer_activate()

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

EngineCombine.mixer_cleanup()

Cleanup the effects of a mixer.

EngineCombine.plot_sweep_stats([axes, …])

Plot sweep_stats to display the convergence with the sweeps.

EngineCombine.plot_update_stats(axes[, …])

Plot update_stats to display the convergence during the sweeps.

EngineCombine.post_update_local(update_data)

Perform post-update actions.

EngineCombine.prepare_svd(theta)

Transform theta into matrix for svd.

EngineCombine.prepare_update()

Prepare self to represent the effective Hamiltonian on sites (i0, i0+1).

EngineCombine.reset_stats()

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

EngineCombine.run()

Run the DMRG simulation to find the ground state.

EngineCombine.set_B(U, S, VH)

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

EngineCombine.sweep([optimize, meas_E_trunc])

One ‘sweep’ of a sweeper algorithm.

EngineCombine.update_LP(U)

Update left part of the environment.

EngineCombine.update_RP(VH)

Update right part of the environment.

EngineCombine.update_local(theta[, …])

Perform bond-update on the sites (i0, i0+1).

Class Attributes and Properties

EngineCombine.DMRG_params

EngineCombine.engine_params

class tenpy.algorithms.dmrg.EngineCombine(psi, model, DMRG_params)[source]

Bases: tenpy.algorithms.dmrg.TwoSiteDMRGEngine

Engine which combines legs into pipes as far as possible.

This engine combines the virtual and physical leg for the left site and right site into pipes. This reduces the overhead of calculating charge combinations in the contractions, but one matvec() is formally more expensive, \(O(2 d^3 \chi^3 D)\).

Deprecated since version 0.5.0: Directly use the TwoSiteDMRGEngine with the DMRG parameter combine=True.

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

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. Only those LP and RP that can be used later should be updated.

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

Return type

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

init_env(model=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.

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

option Sweep.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 Sweep.init_env_data: dict

Dictionary as returned by self.env.get_initialization_data() from get_initialization_data().

option Sweep.orthogonal_to: list of MPSEnvironment

List of other matrix product states to orthogonalize against. Works only for finite systems. This parameter 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.

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.

mixed_svd(theta)[source]

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

The goal is to split theta and truncate it:

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

Without a mixer, this is done by a simple svd and truncation of Schmidt values.

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

Note that the returned S is a general (not diagonal) matrix, with labels 'vL', 'vR'.

Parameters

theta (Array) – The optimized wave function, prepared for svd.

Returns

  • U (Array) – Left-canonical part of theta. Labels '(vL.p0)', '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', '(p1.vR)'.

  • err (TruncationError) – The truncation error introduced.

mixer_activate()[source]

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

option TwoSiteDMRGEngine.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 SingleSiteMixer for the 1-site case.

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

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

  • yaxis (xaxis,) – 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(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_data (dict) – Data computed during the local update, as described in the following list.

  • meas_E_trunc (bool, optional) – Wheter to measure the energy after truncation.

prepare_svd(theta)[source]

Transform theta into matrix for svd.

prepare_update()[source]

Prepare self to represent the effective Hamiltonian on sites (i0, i0+1).

Returns

theta – Current best guess for the ground state, which is to be optimized. Labels 'vL', 'p0', 'vR', 'p1'.

Return type

Array

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

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.

set_B(U, S, VH)[source]

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

Parameters
  • VH (U,) – 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.

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

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

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

Update left part of the environment.

We always update the environment at site i0 + 1: this environment then contains the site where we just performed a local update (when sweeping right).

Parameters

U (Array) – The U as returned by the SVD, with combined legs, labels 'vL.p0', 'vR'.

update_RP(VH)[source]

Update right part of the environment.

We always update the environment at site i0: this environment then contains the site where we just performed a local update (when sweeping left).

Parameters

VH (Array) – The VH as returned by SVD, with combined legs, labels 'vL', '(vR.p1)'.

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

Perform bond-update on the sites (i0, i0+1).

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

  • meas_E_trunc (bool) – Wheter to measure the energy after truncation.

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