hqs_distiller.solver
Module providing LegacyDistiller class.
This module provides the LegacyDistiller class which implements the FEAST algorithm to diagonalize a hermitean, sparse matrix in a user-specified energy window. The FEAST algorithm is based on a rational filter to project trial eigenstates into the subspace spanned by the eigenstates lying in the user-specified energy window. The matrix can then be diagonalized in this subspace using a eigensolver for dense matrices. This implies that the energy window should be chosen small enough to allow for an efficient use of the dense eigensolver. The rational filter is approximated using a Chebyshev expansion of the resolvent matrix, which renders the implementation memory efficient.
Copyright © 2022-2023 HQS Quantum Simulations GmbH. All Rights Reserved.
Class representing solver for partial spectral decompositions.
A solver object can be instantiated using a sparse matrix (scipy.sparse.spmatrix). The instantiated object provides the following user methods:
- get_spectral_decomposition_direct: Solves linear problems using a direct solver.
- get_spectral_decomposition_iterative: Solves linear problems using an iterative solver.
- get_spectral_decomposition_chebyshev: Solves linear problems using a Chebyshev expansion.
A calculation using the Chebyshev expansion depends at the number of trial states and the number of Chebyshev moments used in the projection. The function 'get_number_of_chebyshev_moments' implements a heuristic to set the number of Chebyshev moments used in the projection. It uses the spectral interval for the projection together with the total spectral width of the Hamiltonian and the (estimated) number of states in the spectral interval to determine the number of Chebyshev moments via the formula:
n_chebyshev_moments = C * n_states * spectral_width / spectral_window
NOTE: The constant C is hard-coded to 1/2 at the moment.
The number of trial states is determined by:
n_trial_states = A * n_states
NOTE: The constant A is hard-coded to 3 at the moment.
The user can manually specify the number of Chebyshev moments.
Initialize solver.
After initialization this object provides methods to extract eigenvalue, eigenvector pairs falling into a user-specified energy window.
NOTE: Only for hermitian matrices (real and complex). Not checked in initialization!!!
Arguments:
- matrix (scipy.sparse.spmatrix): Sparse matrix.
- spectral_buffer (float): relative buffer for spectral window. Must be in (0, 0.9). Defaults to 0.05.
- max_n_chebyshev_moments (int): (Maximal) number of Chebyshev moments used. Defaults to 8192.
- verbose (bool): General flag triggering outputs. Defaults to True.
- debug_plots (bool): Flag triggering additional svg outputs. NOTE: verbose and debug_plots must set to True to obtain the svg files. Defaults to False.
Raises:
- ValueError: If matrix is not square matrix.
Get partial spectral decomposition.
NOTE: This method assumes hermiticity of Hamiltonian.
NOTE: Only use for small matrices (linear size ~ 10 k).
FIXME: This method is only for benchmarking purposes.
This method implements the FEAST algorithm using an direct solver (scipy's factorize with SuperLU backend). Memory can be saved by recomputing the LU factors in every iteration, which comes an increased runtime due to the repeated factorization.
Arguments:
- window (tuple[float, float]): Energy window for spectral decomposition.
- states (Optional[np.ndarray]): Initial guess for the trial states. Defaults to None (triggers automatic generation of trial states).
- max_iter (int): Maximum number of iterations for solver. Defaults to 10.
- eps (float): Convergence criterion for solver. Defaults to 1.e-9.
- n_states_limit (int): Maximum allowed number of states in spectral window. Defaults to 512.
- n_initial (int): Number of states used for initial DOS calculation. Defaults to 8.
- verbose (Optional[bool]): Flag triggering output. Defaults to None, which implies that the verbose flag provided during initialization is used.
- n_filter (int): Number of projectors for rational filter. Defaults to 4.
- precompute_factors (bool): Whether to compute the iLU factors once. Defaults to True.
Returns:
tuple[np.ndarray, np.ndarray, np.ndarray]: Tuple of list: Flags indicating convergence of algorithm in each window, list: list of Float arrays containing eigenvalues within each spectral windows. list: list of Complex arrays containing eigenstates within spectral windows.
Get partial spectral decomposition.
NOTE: This method assumes hermiticity of Hamiltonian.
NOTE: Only use for small matrices (linear size ~ 10 k).
NOTE: This method is only for benchmarking purposes.
This method implements the FEAST algorithm using an iterative solver (scipy's gcrotmk) in conjuction with a preconditioner obtained via an incomplete LU (iLU) factorization (scipy's spilu with SuperLU backend). The argument 'fill_factor' determines by how much the iLU factors can increase the non-zero elements of their sparse representation. Memory can be saved by recomputing the iLU's in every iteration, which comes an increased runtime due to the repeated factorization.
Arguments:
- window (tuple[float, float]): Energy window for spectral decomposition.
- states (Optional[np.ndarray]): Initial guess for the trial states. Defaults to None (triggers automatic generation of trial states).
- max_iter (int): Maximum number of iterations for solver. Defaults to 10.
- eps (float): Convergence criterion for solver. Defaults to 1.e-9.
- n_states_limit (int): Maximum allowed number of states in spectral window. Defaults to 512.
- n_initial (int): Number of states used for initial DOS calculation. Defaults to 8.
- verbose (Optional[bool]): Flag triggering output. Defaults to None, which implies that the verbose flag provided during initialization is used.
- n_filter (int): Number of projectors for rational filter. Defaults to 4.
- precompute_factors (bool): Whether to compute the iLU factors once. Defaults to True.
- fill_factor (float): Allowed fill-in for iLU factors. Defaults to 100.
Returns:
tuple[np.ndarray, np.ndarray, np.ndarray]: Tuple of list: Flags indicating convergence of algorithm in each window, list: list of Float arrays containing eigenvalues within each spectral windows. list: list of Complex arrays containing eigenstates within spectral windows.
Get partial spectral decomposition using a chebyshev expansion.
NOTE: This method assumes hermiticity of Hamiltonian.
NOTE: The trial states proposed by the user are augmented with randomly generated trial states to span the subspace to the sub-space dimension auto-detetcted by the algorithm for the given spectral window. The trial states provided by the user are only used if there energy is within the spectral window.
Arguments:
- window (tuple[float, float]): Energy window for spectral decomposition.
- states (Optional[np.ndarray]): Initial guess for eigenstates around spectral window. Defaults to None.
- eps (float): Convergence criterion for solver. Defaults to 1.e-9.
- max_iter (int): Maximum number of iterations for solver. Defaults to 10.
- n_states_limit (int): Maximum allowed number of states in spectral window. Defaults to 64.
- n_initial (int): Number of states used for initial DOS calculation. Defaults to 8.
- verbose (Optional[bool]): Flag triggering output. Defaults to None, which implies that the verbose flag provided during initialization is used.
- n_chebyshev_moments (Optional[int]): Number of Chebyshev moments used for projector. Defaults to None, which determines the number based on provided spectral window.
Returns:
tuple[np.ndarray, np.ndarray, np.ndarray]: Tuple of boolean array indicating convergence for each state in spectral window, float array containing eigenvalues within spectral window. complex array containing eigenstates within spectral window.