LanczosGroundState

  • full name: tenpy.linalg.lanczos.LanczosGroundState

  • parent module: tenpy.linalg.lanczos

  • type: class

class tenpy.linalg.lanczos.LanczosGroundState(H, psi0, params, orthogonal_to=[])[source]

Bases: object

Lanczos algorithm working on npc arrays.

The Lanczos algorithm can finds extremal eigenvalues (in terms of magnitude) along with the corresponding eigenvectors. It assumes that the linear operator H is hermitian. Given a start vector psi0, it generates an orthonormal basis of the Krylov space, in which H is a small tridiagonal matrix, and solves the eigenvalue problem there. Finally, it transform the resulting ground state back into the original space.

Parameters
HNpcLinearOperator-like

A hermitian linear operator. Must implement the method matvec acting on a Array; nothing else required. The result has to have the same legs as the argument.

psi0Array

The starting vector defining the Krylov basis. For finding the ground state, this should be the best guess available. Note that it must not be a 1D “vector”, we are fine with viewing higher-rank tensors as vectors.

paramsdict

Further optional parameters as described in the following table. Add a parameter verbose >=1 to print the used parameters during runtime. The algorithm stops if both criteria for e_tol and p_tol are met or if the maximum number of steps was reached.

key

type

description

N_min

int

Minimum number of steps to perform.

N_max

int

Maximum number of steps to perform.

E_tol

float

Stop if energy difference per step < E_tol

P_tol

float

Tolerance for the error estimate from the Ritz Residual, stop if (RitzRes/gap)**2 < P_tol

min_gap

float

Lower cutoff for the gap estimate used in the P_tol criterion.

N_cache

int

The maximum number of psi to keep in memory during the first iteration. By default, we keep all states (up to N_max). Set this to a number >= 2 if you are short on memory. The penalty is that one needs another Lanczos iteration to determine the ground state in the end, i.e., runtime is large.

reortho

bool

For poorly conditioned matrices, one can quickly loose orthogonality of the generated Krylov basis. If reortho is True, we re-orthogonalize against all the vectors kept in cache to avoid that problem.

cutoff

float

Cutoff to abort if beta (= norm of next vector in Krylov basis before normalizing) is too small. This is necessary if the rank of A is smaller than N_max - then we get a complete basis of the Krylov space, and beta will be zero.

orthogonal_tolist of Array

Vectors (same tensor structure as psi) against which Lanczos will orthogonalize, ensuring that the result is perpendicular to them. (Assumes that the smallest eigenvalue is smaller than 0, which should always be the case if you want to find ground states with Lanczos!)

Notes

I have computed the Ritz residual RitzRes according to http://web.eecs.utk.edu/~dongarra/etemplates/node103.html#estimate_residual. Given the gap, the Ritz residual gives a bound on the error in the wavefunction, err < (RitzRes/gap)**2. The gap is estimated from the full Lanczos spectrum.

Attributes
HNpcLinearOperator-like

The hermitian linear operator.

psi0Array

The starting vector.

orthogonal_tolist of Array

Vectors to orthogonalize against.

N_min, N_max, E_tol, P_tol, N_cache, reortho:

Parameters as described above.

Esndarray, shape(N_max, N_max)

Es[n, :] contains the energies of _T[:n+1, :n+1] in step n.

_Tndarray, shape (N_max + 1, N_max +1)

The tridiagonal matrix representing H in the orthonormalized Krylov basis.

_cutofffloat

See parameter cutoff.

_cachelist of psi0-like vectors

The ONB of the Krylov space generated during the iteration. FIFO (first in first out) cache of at most N_cache vectors.

_result_krylovndarray

Result in the ONB of the Krylov space: ground state of _T.

Methods

run(self)

Find the ground state of H.

run(self)[source]

Find the ground state of H.

Returns
E0float

Ground state energy (estimate).

psi0Array

Ground state vector (estimate).

Nint

Used dimension of the Krylov space, i.e., how many iterations where performed.