Density matrix renormalization group (DMRG)#

DMRG is a variational technique devised to obtain a state that extremizes the expectation value of the Hermitian operator. A typical example is the search for the best MPS approximation of a ground-state of some Hamiltonian in the form of MPO.

A high-level function organizing DMRG simulations is yastn.tn.mps.dmrg_(). See examples at Low energy states of the XX model.

yastn.tn.mps.dmrg_(psi, H, project=None, method='1site', energy_tol=None, Schmidt_tol=None, max_sweeps=1, iterator_step=None, opts_eigs=None, opts_svd=None, precompute=False)[source]#

Perform DMRG sweeps until convergence, starting from MPS psi.

The inner loop sweeps over MPS updating sites from the first site to the last and back, constituting a single iteration. Convergence can be controlled based on energy and/or Schmidt values (which is a more sensitive measure of convergence). The DMRG algorithm sweeps through the lattice at most max_sweeps times or until all convergence measures (with provided tolerance other than the default None) change by less than the provided tolerance during a single sweep.

Outputs iterator if iterator_step is given, which allows inspecting psi outside of dmrg_ function after every iterator_step sweeps.

Parameters:
  • psi (yastn.tn.mps.MpsMpoOBC) – Initial state. It is updated during execution. If psi is not already canonized to the first site, it will be canonized at the start of the algorithm. The output state from dmrg_ is canonized to the first site.

  • H (yastn.tn.mps.MpsMpoOBC | Sequence) – MPO (or a sum of MPOs) to minimize against, see Env().

  • project (Sequence[yastn.tn.mps.MpsMpoOBC | tuple[float, yastn.tn.mps.MpsMpoOBC]]) – Add a penalty to the directions spanned by MPSs in the list. In practice, the penalty is a multiplicative factor that adds a penalty term to the Hamiltonian. As a result, the energy of said MPS rizes and another lowest-energy is targeted by energy minimization. It can be used to find a few low-energy states of the Hamiltonian if the penalty is larger than the energy gap from the ground state. Use [(penalty, MPS), ...] to provide individual penalty for each MPS by hand as a list of tuples (penalty, MPS), where penalty is a number and MPS is an MPS object. If input is a list of MPSs, i.e., [mps, ...], the option uses default penalty=100.

  • method (str | yastn.Method) – DMRG variant to use; options are '1site' or '2site'. Auxlliary class yastn.Method can be used to change the method in between sweeps while the yield gets called after every iterator_step sweeps.

  • energy_tol (float) – Convergence tolerance for the change of energy in a single sweep. By default is None, in which case energy convergence is not checked.

  • Schmidt_tol (float) – Convergence tolerance for the change of Schmidt values on the worst cut/bond in a single sweep. By default is None, in which case Schmidt values convergence is not checked.

  • max_sweeps (int) – Maximal number of sweeps.

  • iterator_step (int) – If int, dmrg_ returns a generator that would yield output after every iterator_step sweeps. The default is None, in which case dmrg_ sweeps are performed immediately.

  • opts_eigs (dict) – options passed to yastn.eigs(). If None, use default {‘hermitian’: True, ‘ncv’: 3, ‘which’: ‘SR’}

  • opts_svd (dict) – Options passed to yastn.linalg.svd_with_truncation() used to truncate virtual spaces (virtual bond dimension of tensors) in method='2site'.

  • precompute (bool) – Controls MPS-MPO-MPS contraction order. If False, use an approach optimal for a single matrix-vector product calculation during iterative Krylov-space building, scaling as O(D^3 M d + D^2 M^2 d^2). If True, uses less optimal contraction order, scaling as O(D^3 M d^2 + D^2 M^2 d^2). However, the latter allows precomputing and reusing part of the diagram during consecutive iterations. Which one is more efficient depends on the parameters. The default is False.

Returns:

  • Generator if iterator_step is not None.

  • DMRG_out(NamedTuple)

    NamedTuple including fields:

    • sweeps number of performed dmrg sweeps.

    • method method used in the last sweep.

    • energy energy after the last sweep.

    • denergy absolut value of energy change in the last sweep.

    • max_dSchmidt norm of Schmidt values change on the worst cut in the last sweep.

    • max_discarded_weight norm of discarded_weights on the worst cut in ‘2site’ procedure.

Single-site DMRG#

In the algorithm we sweep through the MPS, starting from the initial guess psi, optimizing each yastn.Tensor \(A_j\) one by one while keeping all other tensors fixed (alternating least squares). At each step, the best tensor \(A_j\) is found by minimizing the energy of the local effective Hamiltonian \(H_{\rm eff}\). DMRG with method='1site' works with single-site effective Hamiltonians.

# local optimization of the MPS tensor A_j
# using the effective Hamiltonian H_eff
                          ___
                        -|A_j|-
         ___    ___        |        ___    ___    ___
        |___|--|___|-..         ..-|___|--|___|--|___|
          |      |                   |      |      |
         _|_    _|_       _|_       _|_    _|_    _|
H_eff = |___|--|___|-...-|H_j|-...-|___|--|___|--|___|
          |      |         |         |      |      |
         _|_    _|_                 _|_    _|_    _|
        |___|--|___|-..         ..-|___|--|___|--|___|

After the local problem is solved and a new \(A_j\) minimizing \(H_{\rm eff}\) is found, the MPS is orthogonalized to keep it in the canonical form. In the simplest formulation of ‘1site’ DMRG algorithm, the virtual dimension/spaces are fixed.

Two-site DMRG#

The virtual spaces of the MPS can be adjusted by performing DMRG method='2site'. In this approach, we (i) build effective Hamiltonian \(H_{\rm eff}\) of two neighbouring sites, (ii) solve for its ground state, (iii) split the resulting tensor back into two sites with SVD, and then move to next pair. By doing SVD decomposition we can dynamically adjust virtual space of MPS according to parameters in opts_svd.

# local optimization of the MPS tensors A_j and A_j+1 using
# the effective Hamiltonian H_eff

         H_eff of 2-site DMRG
            ___    _____
 ______   -|A_j|--|A_j+1|-    ______
|      |--   |      |      --|      |
|      |                     |      |
|      |    _|_   __|__      |      |
|      |---|H_j|-|H_j+1|-----|      |
|      |     |      |        |      |
|      |                     |      |
|______|--                 --|______|

The next pair will take one of the sites from the previous step while the other one is orthogonalized to maintain canonical form of the MPS.

Projecting out selected MPS#

The optimization can be performed in the restricted subspace, where contributions from some MPSs are projected out. An alternative approach that we utilize is adding penalty terms in the directions of those MPSs. This allows one to search for a few excited states of the Hamiltonian. The list of MPS to project out is given as project=[lower_E_MPS, ...].

Auxilliary mutable Method class#

class yastn.Method(string: str = '')[source]#

Auxiliary mutable method class. It introduces the mechanism to change the method used in yastn.tn.mps.dmrg_(), yastn.tn.mps.tdvp_(), and other generator functions in between consecutive sweeps. Updating the value in place will inject the new value back into the generator.

update_(string)[source]#

Update the method name in place.