Time evolution#

Evolution step#

Time evolution in YASTN operates via the Trotterization of the evolution operator to a set of gates. Each gate acts on a path of sites, single site for a local gate, nearest-neighbor sites for a two-site gate, etc. A single evolution step consists of the application of trotterized gates, together with truncation of the bond dimension after each application of a non-local gate. For a detailed explanation, please refer to the Basic concepts/Time evolution.

The main method for performing a time evolution step is yastn.tn.fpeps.evolution_step_(). This function requires an environment containing the PEPS state (updated in place) and provides control over bond truncation via different ways to calculate bond metric. Classes which support this include:

  • yastn.tn.fpeps.EnvNTU: Employ small local clusters, which can be contracted numerically exactly, resulting in a stable and positively defined bond metric.

  • yastn.tn.fpeps.EnvBP: Employ belief propagation, either to define a bipartite bond metric (equivalent to simple update), or to gauge NTU-like clusters.

  • yastn.tn.fpeps.EnvCTM: Employ CTMRG environment to calculate the bond metric, with information from the entire network, allowing for a full update approach.

  • yastn.tn.fpeps.EnvApproximate: Employ local clusters of sizes beyond exact contraction, utilizing approximate boundary MPS to calculate the bond metric.

yastn.tn.fpeps.evolution_step_(env, gates, opts_svd, method='mpo', fix_metric=0, pinv_cutoffs=(1e-12, 1e-11, 1e-10, 1e-09, 1e-08, 1e-07, 1e-06, 1e-05, 0.0001), max_iter=100, tol_iter=1e-13, initialization='EAT_SVD')[source]#

Perform a single step of PEPS evolution by applying a list of gates. Truncate bond dimension after each application of a two-site gate.

Parameters:
  • env (EnvNTU | EnvCTM | EnvApproximate) – Environment class containing PEPS state (updated in place), and a method to calculate bond metric tensors employed during truncation.

  • gates (yastn.tn.fpeps.Gates) – The gates to be applied to PEPS.

  • opts_svd (dict | Sequence[dict]) – Options passed to yastn.linalg.svd_with_truncation() which are used to initially truncate the bond dimension before it is further optimized iteratively. In particular, it fixes bond dimensions (potentially, sectorial). It is possible to provide a list of dicts (with decreasing bond dimensions), in which case the truncation is done gradually in a few steps.

  • method (str) – For 'NN', split multi-site gates into a series of nearest-neighbor gates. Otherwise, apply mpo-gate first, and then sequentially truncate enlarged bonds (the default).

  • fix_metric (int | None) – Error measure of the metric tensor is a sum of: the norm of its non-hermitian part and absolute value of the most negative eigenvalue (if present). If fix_metric is a number, replace eigenvalues smaller than the error_measure with fix_metric * error_measure. Sensible values of fix_metric are \(0\) and \(1\). The default is \(0\). If None, do not perform eigh to test for negative eigenvalues and do not fix metric.

  • pinv_cutoffs (Sequence[float] | float) – List of pseudo-inverse cutoffs. The one that gives the smallest truncation error is used during iterative optimizations and EAT initialization.

  • max_iter (int) – The maximal number of iterative steps for each truncation optimization.

  • tol_iter (int) – Tolerance of truncation_error to stop iterative optimization.

  • initialization (str) – Tested initializations of iterative optimization. The one resulting in the smallest error is selected. Possible options are ‘SVD’ (svd initialization only), ‘EAT’ (EAT optimization only), ‘SVD_EAT’ (tries both).

Returns:

List of Namedtuple containing fields:
  • bond bond where the gate is applied.

  • truncation_error relative norm of the difference between untruncated and truncated bonds, calculated in metric specified by env.

  • best_method initialization/optimization method giving the best truncation_error. Possible values are ‘eat’, ‘eat_opt’, ‘svd’, ‘svd_opt’.

  • nonhermitian_part norm of the non-hermitian part of the bond metric, normalized by the bond metric norm. Estimator of metric error.

  • min_eigenvalue the smallest bond metric eigenvalue, normalized by the bond metric norm. Can be negative which gives an estimate of metric error.

  • wrong_eigenvalues a fraction of bond metrics eigenvalues that were below the error threshold; and were modified according to fix_metric argument.

  • eat_metric_error error of approximating metric tensor by a product via SVD-1 in EAT initialization.

  • truncation_errors dict with truncation_error-s for all tested initializations/optimizations.

  • iterations dict with number of iterations to converge iterative optimization.

  • pinv_cutoffs dict with optimal pinv_cutoffs for methods where pinv appears.

Return type:

list[Evolution_out(NamedTuple)]

Auxiliary functions, such as yastn.tn.fpeps.accumulated_truncation_error(), assist with tracking cumulative truncation error across multiple evolution steps, facilitating error analysis in time evolution simulations.

yastn.tn.fpeps.accumulated_truncation_error(infoss, statistics='mean')[source]#

Return accumulated truncation error \(\Delta\) calcuated from evolution output statistics.

\(\Delta = \sum_{steps} statistics_{bond} [\sum_{gate \in bond} truncation\_error(gate, step)]\), where statistics is mean or max.

Gives an estimate of errors accumulated during time evolution.

Parameters:
  • infoss (Sequence[Sequence[Evolution_out]]) – list of outputs of evolution_step_().

  • statistics (str) – ‘max’ or ‘mean’, whether to take the maximal value or a mean over the bonds in the lattice.

Example

infoss = []
for step in range(number_of_steps):
    infos = fpeps.evolution_step_(env, gates, opt_svd)
    infoss.append(infos)

# Accumulated runcation error
Delta = fpeps.accumulated_truncation_error(infoss)

To apply a gate to a PEPS, without truncating the bond dimension, one can use yastn.tn.fpeps.Peps.apply_gate_(). Truncation of the bond dimension can be done later, e.g., via yastn.tn.fpeps.truncate_(). Both operations modify the PEPS in place.

Peps.apply_gate_(gate)[source]#

Apply gate to PEPS state psi, changing respective tensors in place.

fpeps.truncate_(opts_svd, bond=None, fix_metric=0, pinv_cutoffs=(1e-12, 1e-11, 1e-10, 1e-09, 1e-08, 1e-07, 1e-06, 1e-05, 0.0001), max_iter=100, tol_iter=1e-13, initialization='EAT_SVD')#

Truncate virtual bond dimensions of PEPS.

Parameters:
  • env (EnvNTU | EnvCTM | EnvApproximate) – Environment class containing PEPS state (updated in place), and a method to calculate bond metric tensors employed during truncation.

  • opts_svd (dict | Sequence[dict]) – Options passed to yastn.linalg.svd_with_truncation() which are used to initially truncate the bond dimension before it is further optimized iteratively. In particular, it fixes bond dimensions (potentially, sectorial). It is possible to provide a list of dicts (with decreasing bond dimensions), in which case the truncation is done gradually in a few steps.

  • bond (tuple[Site, Site] | None) – Specify single bond to be truncated. If None, truncate all bonds. The default is None.

  • fix_metric (int | None) – Error measure of the metric tensor is a sum of: the norm of its non-hermitian part and absolute value of the most negative eigenvalue (if present). If fix_metric is a number, replace eigenvalues smaller than the error_measure with fix_metric * error_measure. Sensible values of fix_metric are \(0\) and \(1\). The default is \(0\). If None, do not perform eigh to test for negative eigenvalues and do not fix metric.

  • pinv_cutoffs (Sequence[float] | float) – List of pseudo-inverse cutoffs. The one that gives the smallest truncation error is used during iterative optimizations and EAT initialization.

  • max_iter (int) – The maximal number of iterative steps for each truncation optimization.

  • tol_iter (int) – Tolerance of truncation_error to stop iterative optimization.

  • initialization (str) – Tested initializations of iterative optimization. The one resulting in the smallest error is selected. Possible options are ‘SVD’ (svd initialization only), ‘EAT’ (EAT optimization only), ‘SVD_EAT’ (tries both).

Returns:

Namedtuple containing fields:
  • bond bond where the gate is applied.

  • truncation_error relative norm of the difference between untruncated and truncated bonds, calculated in metric specified by env.

  • best_method initialization/optimization method giving the best truncation_error. Possible values are ‘eat’, ‘eat_opt’, ‘svd’, ‘svd_opt’.

  • nonhermitian_part norm of the non-hermitian part of the bond metric, normalized by the bond metric norm. Estimator of metric error.

  • min_eigenvalue the smallest bond metric eigenvalue, normalized by the bond metric norm. Can be negative which gives an estimate of metric error.

  • wrong_eigenvalues a fraction of bond metrics eigenvalues that were below the error threshold; and were modified according to fix_metric argument.

  • eat_metric_error error of approximating metric tensor by a product via SVD-1 in EAT initialization.

  • truncation_errors dict with truncation_error-s for all tested initializations/optimizations.

  • iterations dict with number of iterations to converge iterative optimization.

  • pinv_cutoffs dict with optimal pinv_cutoffs for methods where pinv appears.

Return type:

Evolution_out(NamedTuple)

Gates#

yastn.tn.fpeps.evolution_step_() takes a list of gates to be applied on PEPS

class yastn.tn.fpeps.Gate(G: tuple = None, sites: tuple = None)[source]#

Gate to be applied on Peps state.

G contains operators for respective sites.

Operator can be given in the form of an MPO (yastn.tn.mps.MpsMpoOBC) of the same length as the number of provided sites. Sites should form a continuous path in the two-dimensional PEPS lattice. The fermionic order of MPO should be linear, with the first MPO site being first in the fermionic order, irrespective of the provided sites.

For a two-site operator acting on sites beyond nearest neighbor, it can be provided as G with two elements, and sites forming a path between the end sites where the provided elements of G will act.

It is also possible to provide G as a list of tensots. In this case, the convention of legs is (ket, bra, virtual_0, virtual_1) – i.e., the first two legs are always physical (operator) legs. For one site, there are no virtual legs. For two or more sites, the first and last elements of G have one virtual leg (3 in total). For three sites or more, the middle elements of G have two virtual legs connecting, respectively, to the preceding and following elements of G.

Create new instance of Gate(G, sites)

An auxiliary function yastn.tn.fpeps.gates.distribute() distribute a set of gates homogeneously over the entire lattice. By default, it complements gates with their adjoint (gates repeated in reverse order), forming a 2nd order approximation for a small time step – this assumes that individual gates take half of the time step.

yastn.tn.fpeps.gates.distribute(geometry, gates_nn=None, gates_local=None, symmetrize=True) Sequence[Gate][source]#

Distributes gates homogeneous over the lattice.

Parameters:
  • geometry (yastn.tn.fpeps.SquareLattice | yastn.tn.fpeps.CheckerboardLattice | yast.tn.fpeps.Peps) – Geometry of PEPS lattice. Can be any structure that includes geometric information about the lattice, like the Peps class.

  • gates_nn (Gate | Sequence[Gate]) – Nearest-neighbor gate, or a list of gates, to be distributed over all unique lattice bonds.

  • gates_local (Gate | Sequence[Gate]) – Local gate, or a list of local gates, to be distributed over all unique lattice sites.

  • symmetrize (bool) – Whether to iterate through provided gates forward and then backward, resulting in a 2nd order method. In that case, each gate should correspond to half of the desired timestep. The default is True.

Predefined gates#

Some predefined gates can be found in yastn.tn.fpeps.gates, including

yastn.tn.fpeps.gates.gate_nn_hopping(t, step, I, c, cdag, bond=None) Gate[source]#

Nearest-neighbor gate \(G = \exp(-step \cdot H)\) for \(H = -t \cdot (c^\dagger_1 c_2 + c^\dagger_2 c_1)\)

\(G = I + (\cosh(x) - 1) (n_1 h_2 + h_1 n_2) + \sinh(x) (c^\dagger_1 c_2 + c^\dagger_2 c_1)\), where \(x = t \cdot step\)

yastn.tn.fpeps.gates.gate_nn_Ising(J, step, I, X, bond=None) Gate[source]#

Nearest-neighbor gate \(G = \exp(-step \cdot H)\) for \(H = J X_1 X_2\), where \(X\) is a Pauli matrix.

\(G = \cosh(x) I - \sinh(x) X_1 X_2\), where \(x = step \cdot J\).

yastn.tn.fpeps.gates.gate_nn_Heisenberg(J, step, I, Sz, Sp, Sm, bond=None) Gate[source]#

Nearest-neighbor gate \(G = \exp(-step \cdot H)\) for \(H = J (S_z S_z + 2 * S_+ S_- + 2 * S_- S_+)\), for spin operators \(S_z, S_+, S_-\).

yastn.tn.fpeps.gates.gate_local_Coulomb(mu_up, mu_dn, U, step, I, n_up, n_dn, site=None) Gate[source]#

Local gate \(\exp(-step \cdot H)\) for \(H = U \cdot (n_{up} - I / 2) \cdot (n_{dn} - I / 2) - \mu_{up} \cdot n_{up} - \mu_{dn} \cdot n_{dn}\)

We ignore a constant \(U / 4\) in the above Hamiltonian.

yastn.tn.fpeps.gates.gate_local_occupation(mu, step, I, n, site=None) Gate[source]#

Local gate \(G = \exp(-step \cdot H)\) for \(H = -\mu \cdot n\)

\(G = I + n \cdot (\exp(\mu \cdot step) - 1)\)

yastn.tn.fpeps.gates.gate_local_field(h, step, I, X, site=None) Gate[source]#

Local gate \(G = \exp(-step \cdot H)\) for \(H = -h \cdot X\) where \(X\) is a Pauli matrix.

\(G = \cosh(h \cdot step) \cdot I + \sinh(h \cdot step) \cdot X\)

yastn.tn.fpeps.gates.gate_nn_tJ(J, tu, td, muu0, muu1, mud0, mud1, step, I, cu, cpu, cd, cpd, bond=None) Gate[source]#

Nearest-neighbor gate \(G = \exp(-step \cdot H)\) for \(H = -t \sum_{\sigma} (c_{0,\sigma}^\dagger c_{1,\sigma} + c_{1,\sigma}^\dagger c_{0,\sigma}) + J (S_i \cdot S_j - \frac{n_i n_j}{4}) - \sum_{i, \sigma} \mu_{i,\sigma} n_{i,\sigma}\)

Custom gates#

Auxiliary functions help generate gates by exponentiating hermitian Hamiltonians, generate two-site Hamiltonian ensuring proper fermionic order, etc.

yastn.tn.fpeps.fkron(A, B, sites=(0, 1), merge=True)[source]#

Fermionic kron; auxiliary function for gate generation. Returns a Kronecker product of two local operators, A and B, including swap-gate (fermionic string) to handle fermionic operators.

Calculate A_0 B_1 for sites==(0, 1), and A_1 B_0 for sites==(1, 0), i.e., operator B acts first (relevant when operators are fermionically non-trivial).

Site 0 is assumed to be first in fermionic orders.

If merge, returns equivalent of ncon([A, B], [(-0, -1), (-2, -3)]) for sites==(0, 1), and ncon([A, B], [(-2, -3), (-0, -1)]) for sites==(1, 0), with proper operator order and swap gate applied.:

   1     3
   |     |
┌──┴─────┴──┐
|           |
└──┬─────┬──┘
   |     |
   0     2
yastn.tn.fpeps.gates.gate_nn_exp(step, I, H, bond=None) Gate[source]#

Gate exp(-step * H) for a hermitian Hamiltonian H, consistent with leg order of yastn.tn.fpeps.gates.fkron(). Add 0 * I to the Hamiltonian to avoid situation, where some blocks are missing in the Hamiltonian.

yastn.tn.fpeps.gates.gate_local_exp(step, I, H, site=None) Gate[source]#

Gate exp(-step * H) for local Hamiltonian H. Add 0 * I to the Hamiltonian to avoid situation, where some blocks are missing in the Hamiltonian.

yastn.tn.fpeps.gates.decompose_nn_gate(Gnn, bond=None) Gate[source]#

Auxiliary function to generate Gate by cutting a two-site operator, using SVD, into two local operators with the connecting legs.