SingleSiteMixer

Inheritance Diagram

Inheritance diagram of tenpy.algorithms.dmrg.SingleSiteMixer

Methods

SingleSiteMixer.__init__(options)

Initialize self.

SingleSiteMixer.perturb_svd(engine, theta, …)

Mix extra terms to theta and perform an SVD.

SingleSiteMixer.subspace_expand(engine, …)

Expand the MPS subspace, to allow the bond dimension to increase.

SingleSiteMixer.update_amplitude(sweeps)

Update the amplitude, possibly disable the mixer.

class tenpy.algorithms.dmrg.SingleSiteMixer(options)[source]

Bases: tenpy.algorithms.dmrg.Mixer

Mixer for single-site DMRG.

Performs a subspace expansion following [Hubig2015].

perturb_svd(engine, theta, i0, move_right, next_B)[source]

Mix extra terms to theta and perform an SVD.

We calculate the left and right reduced density matrix using the mixer (which might include applications of H). These density matrices are diagonalized and truncated such that we effectively perform a svd for the case mixer.amplitude=0.

Parameters
  • engine (Engine) – The DMRG engine calling the mixer.

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

  • i0 (int) – The site index where theta lives.

  • move_right (bool) – Whether we move to the right (True) or left (False).

  • next_B (Array) – The subspace expansion requires to change the tensor on the next site as well. If move_right, it should correspond to engine.psi.get_B(i0+1, form='B'). If not move_right, it should correspond to engine.psi.get_B(i0-1, form='A').

Returns

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

  • S (1D ndarray) – (Perturbed) singular values on the new bond (between theta and next_B).

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

  • err (TruncationError) – The truncation error introduced.

subspace_expand(engine, theta, i0, move_right, next_B)[source]

Expand the MPS subspace, to allow the bond dimension to increase.

This is the subspace expansion following [Hubig2015].

Parameters
  • engine (SingleSiteDMRGEngine | TwoSiteDMRGEngine) – ‘Engine’ for the DMRG algorithm

  • theta (Array) – Optimized guess for the ground state of the effective local Hamiltonian.

  • i0 (int) – Site index at which the local update has taken place.

  • move_right (bool) – Whether the next i0 of the sweep will be right or left of the current one.

  • next_B (Array) – The subspace expansion requires to change the tensor on the next site as well. If move_right, it should correspond to engine.psi.get_B(i0+1, form='B'). If not move_right, it should correspond to engine.psi.get_B(i0-1, form='A').

Returns

  • theta – Local MPS tensor at site i0 after subspace expansion.

  • next_B – MPS tensor at site i0+1 or i0-1 (depending on sweep direction) after subspace expansion.

update_amplitude(sweeps)[source]

Update the amplitude, possibly disable the mixer.

Parameters

sweeps (int) – The number of performed sweeps, to check if we need to disable the mixer.

Returns

mixer – Returns self if we should continue mixing, or None, if the mixer should be disabled.

Return type

Mixer | None