Time evolution#

Evolution step#

Time evolution in YASTN operates via the Trotterization of evolution operator to a set of nearest-neighbor and local gates. A single evolution step consists of the application of trotterized gates, together with truncation of the bond dimension after each application of a nearest-neighbor gate. For a detailed explanation, please refer to the Time Evolution section .

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

  • :class:`yastn.tn.fpeps.EnvNTU`: Suitable for local truncation with faster simulation.

  • :class:`yastn.tn.fpeps.EnvCTM`: Used for fast full update schemes with high accuracy.

  • :class:`yastn.tn.fpeps.EnvApproximate`: A lower-cost alternative for approximate updates on large clusters (Marek ?).

yastn.tn.fpeps.evolution_step_(env, gates, opts_svd, symmetrize=True, 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-15, 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 | …) – Environment class containing PEPS state (updated in place), and a method to calculate bond metric tensors employed during truncation.

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

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

  • 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 not None, replace eigenvalues smaller than the error_measure by 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)

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(sum_bond truncation_error(bond, step))} 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 of a mean over 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)

Gates#

Gates classes are organized as

class yastn.tn.fpeps.Gates(nn: list = (), local: list = ())[source]#

List of nearest-neighbor and local operators to be applied to PEPS during evolution_step.

Create new instance of Gates(nn, local)

class yastn.tn.fpeps.Gate_nn(G0: tuple | None = None, G1: tuple | None = None, bond: tuple | None = None)[source]#

G0 should be before G1 in the fermionic and lattice orders. The third legs of G0 and G1 are auxiliary legs connecting them into a two-site operator.

If a bond is None, this is a general operator. Otherwise, the bond can carry information where it should be applied (potentially, after fixing the order mismatches).

Create new instance of Gate_nn(G0, G1, bond)

class yastn.tn.fpeps.Gate_local(G: tuple | None = None, site: tuple | None = None)[source]#

G is a local operator with ndim==2.

If site is None, this is a general operator. Otherwise, the site can carry information where it should be applied.

Create new instance of Gate_local(G, site)

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_nn[source]#

Nearest-neighbor gate G = exp(-step * H) for H = -t * (cdag_1 c_2 + cdag_2 c_1)

G = I + (cosh(x) - 1) * (n_1 h_2 + h_1 n_2) + sinh(x) * (cdag_1 c_2 + cdag_2 c_1)

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

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

G = cosh(x) I - sinh(x) X_1 X_2; x = step * J

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

Local gate exp(-step * H) for H = U * (n_up - I / 2) * (n_dn - I / 2) - mu_up * n_up - mu_dn * 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_local[source]#

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

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

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

Local gate G = exp(-step * H) for H = -h * X

G = cosh(h * step) * I + np.sinh(h * step) * X

An auxiliary function yastn.tn.fpeps.gates.distribute() distribute a set of gates homogeneously over the entire lattice.

yastn.tn.fpeps.gates.distribute(geometry, gates_nn=None, gates_local=None) Gates[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.

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

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