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_()
. For examples, see DMRG.
- 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)[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 None) change by less than the provided tolerance during a single sweep.Outputs iterator if
iterator_step
is given, which allows inspectingpsi
outside ofdmrg_
function after everyiterator_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)
, wherepenalty
is a number andMPS
is an MPS bject. If input is a list of MPSs, i.e.,[mps, ...]
, the option uses defaultpenalty=100
.method (str) – DMRG variant to use; options are ‘1site’ or ‘2site’.
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. Default is None, in which casedmrg_
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) inmethod='2site'
.
- Returns:
Generator if iterator_step is not None.
DMRG_out(NamedTuple) –
NamedTuple including fields:
sweeps
number of performed dmrg sweeps.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. This can be useful when searching for
excited states. List of MPS to project out is given as project=[lower_E_MPS, ...]
.