quspin.tools.lanczos.FTLM_static_iteration

quspin.tools.lanczos.FTLM_static_iteration(O_dict, E, V, Q_T, beta=0)[source]

Calculate iteration for Finite-Temperature Lanczos method.

Here we give a brief overview of this method based on arXiv:1111.5931.

One would naively think that it would require full diagonalization to calculate thermodynamic expectation values for a quantum system as one has to fully diagonalize the Hamiltonian to evaluate:

\[\langle O\rangle_\beta = \frac{1}{Z}Tr\left(e^{-\beta H}O\right)\]

with the partition function defined as: \(Z=Tr\left(e^{-\beta H}\right)\). The idea behind the Finite-Temperature Lanczos Method (FTLM) is to use quantum typicality as well as Krylov subspaces to simplify this calculation. Typicality states that the trace of an operator can be approximated as an average of that same operator with random vectors in the Hilbert-space sampled with the Harr measure. As a corollary, it is known that the fluctuations of this average for any finite sample set will converge to 0 as the size of the Hilbert space increases. Mathematically this is expressed as:

\[\frac{1}{\dim\mathcal{H}}Tr\left(e^{-\beta H}O\right)\approx \frac{1}{N_r}\sum_r\langle r| e^{-\beta H}O |r\rangle\]

where \(|r\rangle\) is a random state from the Harr measure of hilbert space \(\mathcal{H}\) if the Hamiltonian. An issue can occur when the temperature goes to zero as the overlap \(\langle r| e^{-\beta H}O |r\rangle\) will be quite small for most states \(|r\rangle\). Hence, this will require more random realizations to converge. Therefore, this method is really only useful at high temperatures. Next, the idea is to use Lanczos to approximate the matrix exponential. The eigenstates from the lanczos basis can effectively be inserted as an identity operator:

\[\frac{1}{\dim\mathcal{H}}Tr\left(e^{-\beta H}O\right)\approx \frac{1}{N_r}\sum_r\langle r|e^{-\beta H}O|r\rangle\approx \frac{1}{N_r}\sum_r\sum_{i=1}^m e^{-\beta\epsilon^{(r)}_i}\langle r|\psi^{(r)}_i\rangle\langle\psi^{(r)}_i|O|r\rangle = \frac{1}{N_r}\sum_r \langle O\rangle_r \equiv \overline{\langle O\rangle_r}\]

Now going back to the thermal expecation value, we can use the expression above to calculate \(\frac{1}{Z}Tr\left(e^{-\beta H}O\right)\) by noting that the partition function is simply the expecation value of the identity operator: \(Z=Tr\left(e^{-\beta H}I\right)\) and hence the thermal expecation value is approximated by:

\[\langle O\rangle_\beta \approx \frac{\overline{\langle O\rangle_r}}{\overline{\langle I\rangle_r}}\]

The idea behind this function is to generate the the expecation value \(\langle O\rangle_r\) and \(\langle I\rangle_r\) for a lanczos basis generated from an initial state \(|r\rangle\). Therefore if the user would like to calculate the thermal expecation value all one has to do call this function for each lanczos basis generated from a random state \(|r\rangle\).

Parameters:
O_dictdictionary of Python Objects

These Objects must have a ‘dot’ method that calculates a matrix vector product on a numpy.ndarray[:], the effective shape of these objects should be (n,n).

Earray_like, (m,)

Eigenvalues for the Krylov projection of some operator.

Varray_like, (m,m)

Eigenvectors for the Krylov projection of some operator.

Q_Titerator over rows of Q_T

generator or ndarray that contains the lanczos basis associated with E, and V.

betascalar/array_like, any shape

Inverse temperature values to evaluate.

Returns:
Result_dict: dictionary

A dictionary storying the results for a single iteration of the FTLM. The results are stored in numpy.ndarrays that have the same shape as beta. The keys of Result_dict are the same as the keys in O_dict and the values associated with the given key in Result_dict are the expectation values for the operator in O_dict with the same key.

I_expt: numpy.ndarray, same shape as beta

The expecation value of the identity operator for each beta.

Notes

  • The amount of memory used by this function scales like: \(nN_{op}\) with \(n\) being the size of the full Hilbert space and \(N_{op}\) is the number of input operators.

  • FTLM does not converge very well at low temperatures, see function for low-temperature lanczos iterations.

  • One has to be careful as typicality only applies to the trace operation over the entire Hilbert space. Using symmetries is possible, however it requires the user to keep track of the weights in the different sectors.

Examples

>>> beta = numpy.linspace(0,10,101)
>>> E, V, Q_T = lanczos_full(H,v0,20)
>>> Res,Id = FTLM_static_iteration(Obs_dict,E,V,Q_T,beta=beta)