Spectra

The methods in this Spectra module uses our another package eMaTe, which is available at https://github.com/stdogpkg/emate .

eMaTe is a python package implemented in tensorflow which the main goal is provide useful methods capable of estimate spectral densities and trace functions of large sparse matrices.

Trace Functions

The core module responsible to calc trace functions.

Given a semi-positive definite matrix \(A \in \mathbb R^{|V|\times|V|}\), which has the set of eigenvalues given by \(\{\lambda_i\}\) a trace of a matrix function is given by

\[\mathrm{tr}(f(A)) = \sum\limits_{i=0}^{|V|} f(\lambda_i)\]

The methods for calculating such traces functions have a cubic computational complexity lower bound, \(O(|V|^3)\). Therefore, it is not feasible for  large networks. One way to overcome such computational complexity it is use stochastic approximations combined with a mryiad of another methods to get the results with enough accuracy and with a small computational cost. The methods available in this module uses the Sthocastic Lanczos Quadrature, a procedure proposed in the work made by Ubaru, S. et.al. [1] (you need to cite them).

References

[1] Ubaru, S., Chen, J., & Saad, Y. (2017). Fast Estimation of tr(f(A)) via Stochastic Lanczos Quadrature. SIAM Journal on Matrix Analysis and Applications, 38(4), 1075-1099.

[2] Hutchinson, M. F. (1990). A stochastic estimator of the trace of the influence matrix for laplacian smoothing splines. Communications in Statistics-Simulation and Computation, 19(2), 433-450.

stdog.spectra.trace_function.slq(A, num_vecs, num_steps, trace_function, device='/gpu:0', precision=32, random_factory=<function radamacher>, parallel_iterations=10, swap_memory=False, infer_shape=False)

Compute the approxiamted value of a given trace function using the sthocastic Lanczos quadrature using Radamacher’s random vectors.

Parameters:
  • A (scipy sparse matrix) – The semi-positive definite matrix
  • num_vecs (int) – Number of random vectors in oder to aproximate the trace
  • num_steps (int) – Number of Lanczos steps
  • trace_function (function) –

    A function like

    def trace_function(eig_vals)
        *tensorflow ops
        return result
    
  • precision (int) –
    1. Single or (64) double precision
Returns:

  • f_estimation (float) – The approximated value of the given trace function
  • gammas (array of floats) – See [1] for more details

References

[1] Ubaru, S., Chen, J., & Saad, Y. (2017). Fast Estimation of tr(f(A)) via Stochastic Lanczos Quadrature. SIAM Journal on Matrix Analysis and Applications, 38(4), 1075-1099.

stdog.spectra.trace_function.entropy(L_sparse, num_vecs=100, num_steps=50, device='/gpu:0')[source]

Compute the spectral entropy

\[ \begin{align}\begin{aligned}\sum\limits_{i=0}^{|V|} f(\lambda_i)\\f(\lambda) = \begin{cases} -\lambda \log_2\lambda \ \ if \ \lambda > 0; \newline 0,\ \ otherwise \end{cases}\end{aligned}\end{align} \]
Parameters:
  • L (sparse matrix) –
  • num_vecs (int) – Number of random vectors used to approximate the trace using the Hutchison’s trick [1]
  • num_steps (int) – Number of Lanczos steps or Chebyschev’s moments
  • device (str) – “/cpu:int” our “/gpu:int”
Returns:

approximated_spectral_entropy

Return type:

float

References

1 - Ubaru, S., Chen, J., & Saad, Y. (2017). Fast Estimation of tr(f(A)) via Stochastic Lanczos Quadrature. SIAM Journal on Matrix Analysis and Applications, 38(4), 1075-1099.

2 - Hutchinson, M. F. (1990). A stochastic estimator of the trace of the influence matrix for laplacian smoothing splines. Communications in Statistics-Simulation and Computation, 19(2), 433-450.

stdog.spectra.trace_function.estrada_index(L, num_vecs=100, num_steps=50, device='/gpu:0')[source]

Given the Laplacian matrix \(L \in \mathbb R^{|V|\times |V|}\) s.t. for all \(v \in \mathbb R^{|V|}\) we have \(v^T L v > 0\) the Estrada Index is given by

\[\mathrm{tr}\exp(L) = \sum\limits_{i=0}^{|V|} e^{\lambda_i}\]
Parameters:
  • L (sparse matrix) –
  • num_vecs (int) – Number of random vectors used to approximate the trace using the Hutchison’s trick [1]
  • num_steps (int) – Number of Lanczos steps or Chebyschev’s moments
  • device (str) – “/cpu:int” our “/gpu:int”
Returns:

approximated_estrada_index

Return type:

float

References

1 - Ubaru, S., Chen, J., & Saad, Y. (2017). Fast Estimation of tr(f(A)) via Stochastic Lanczos Quadrature. SIAM Journal on Matrix Analysis and Applications, 38(4), 1075-1099.

2 - Hutchinson, M. F. (1990). A stochastic estimator of the trace of the influence matrix for laplacian smoothing splines. Communications in Statistics-Simulation and Computation, 19(2), 433-450.

Spectral Density

The Kernel Polynomial Method can estimate the spectral density of large sparse Hermitan matrices with a computational cost almost linear. This method combines three key ingredients: the Chebyshev expansion + the stochastic trace estimator + kernel smoothing.

stdog.spectra.dos.kpm(H, num_moments, num_vecs, extra_points, precision=32, lmin=None, lmax=None, epsilon=0.01, device='/gpu:0', swap_memory_while=False)

Kernel Polynomial Method using a Jackson’s kernel.

Parameters:
  • H (scipy sparse matrix) – The Hermitian matrix
  • num_moments (int) –
  • num_vecs (int) – Number of random vectors in oder to aproximate the trace
  • extra_points (int) –
  • precision (int) – Single or double precision
  • limin (float, optional) – The smallest eigenvalue
  • lmax (float) – The highest eigenvalue
  • epsilon (float) – Used to rescale the matrix eigenvalues into the interval [-1, 1]
Returns:

  • ek (array of floats) – An array with num_moments + extra_points approximated “eigenvalues”
  • rho (array of floats) – An array containing the densities of each “eigenvalue”

References

[1] Wang, L.W., 1994. Calculating the density of states and optical-absorption spectra of large quantum systems by the plane-wave moments method. Physical Review B, 49(15), p.10154.

[2] Hutchinson, M.F., 1990. A stochastic estimator of the trace of the influence matrix for laplacian smoothing splines. Communications in Statistics-Simulation and Computation, 19(2), pp.433-450.