solver

Eigensolvers with a few extra computation methods

The Solver class is the main interface for dealing with eigenvalue problems. It is made to work specifically with pybinding’s Model objects, but it may use any eigensolver algorithm under the hood.

A few different algorithms are provided out of the box: the lapack(), arpack() and feast() functions return concrete Solver implementation using the LAPACK, ARPACK and FEAST algorithms, respectively.

The Solver may easily be extended with new eigensolver algorithms. All that is required is a function which takes a Hamiltonian matrix and returns the computed eigenvalues and eigenvectors. See _SolverPythonImpl for example.

Classes

Solver(impl)

Computes the eigenvalues and eigenvectors of a Hamiltonian matrix

_SolverPythonImpl(solve_func, model, **kwargs)

Python eigensolver implementation

Functions

arpack(model, k[, sigma])

ARPACK Solver implementation for sparse matrices

dacp(model[, window, random_vectors, ...])

pyDACP Solver implementation for DACP method matrices

feast(model, energy_range, initial_size_guess)

FEAST Solver implementation for sparse matrices

lapack(model, **kwargs)

LAPACK Solver implementation for dense matrices

class Solver(impl: Solver | _SolverPythonImpl)

Computes the eigenvalues and eigenvectors of a Hamiltonian matrix

This the common interface for various eigensolver implementations. It should not be created directly, but via the specific functions: lapack(), arpack() and feast(). Those functions will set up their specific solver strategy and return a properly configured Solver object.

calc_bands(k0: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], k1: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], *ks: Iterable[_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]], step: float = 0.1, point_labels: List[str] | None = None) Bands

Calculate the band structure on a path in reciprocal space

Parameters:
k0, k1, *ksarray_like

Points in reciprocal space which form the path for the band calculation. At least two points are required.

stepfloat, optional

Calculation step length in reciprocal space units. Lower step values will return more detailed results.

point_labelsList[str], optional

The point_labels for plots

Returns:
Bands
calc_bands_area(k_area: Area) BandsArea

Calculate the band structure on a path in reciprocal space

Parameters:
k_area~pybinding.Area

Points in reciprocal space which form the Area for the band calculation. At least two points are required.

Returns:
BandsArea
calc_bands_path(k_path: Path) Bands

Calculate the band structure on a path in reciprocal space

Parameters:
k_path~pybinding.Path

Points in reciprocal space which form the path for the band calculation. At least two points are required.

Returns:
Bands
calc_dos(energies: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], broadening: float) Series

Calculate the density of states as a function of energy

\[\text{DOS}(E) = \frac{1}{c \sqrt{2\pi}} \sum_n{e^{-\frac{(E_n - E)^2}{2 c^2}}}\]

for each \(E\) in energies, where \(c\) is broadening and \(E_n\) is eigenvalues[n].

Parameters:
energiesarray_like

Values for which the DOS is calculated.

broadeningfloat

Controls the width of the Gaussian broadening applied to the DOS.

Returns:
Series
calc_eigenvalues(map_probability_at: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes] | None = None) Eigenvalues

Return an Eigenvalues result object with an optional probability colormap

While the eigenvalues property returns the raw values array, this method returns a result object with more data. In addition to the energy states, this result may show a colormap of the probability density for each state at a single position.

Parameters:
map_probability_atarray_like, optional

Cartesian position where the probability density of each energy state should be calculated.

Returns:
Eigenvalues
calc_ldos(energies: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], broadening: float, position: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], sublattice: str = '', reduce: bool = True) Series

Calculate the local density of states as a function of energy at the given position

\[\text{LDOS}(E) = \frac{1}{c \sqrt{2\pi}} \sum_n{|\Psi_n(r)|^2 e^{-\frac{(E_n - E)^2}{2 c^2}}}\]

for each \(E\) in energies, where \(c\) is broadening, \(E_n\) is eigenvalues[n] and \(r\) is a single site position determined by the arguments position and sublattice.

Parameters:
energiesarray_like

Values for which the DOS is calculated.

broadeningfloat

Controls the width of the Gaussian broadening applied to the DOS.

positionarray_like

Cartesian position of the lattice site for which the LDOS is calculated. Doesn’t need to be exact: the method will find the actual site which is closest to the given position.

sublatticestr

Only look for sites of a specific sublattice, closest to position. The default value considers any sublattice.

reducebool

This option is only relevant for multi-orbital models. If true, the resulting LDOS will summed over all the orbitals at the target site and the result will be a 1D array. If false, the individual orbital results will be preserved and the result will be a 2D array with shape == (energy.size, num_orbitals).

Returns:
Series
calc_probability(n: int | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], reduce: float = 1e-05) StructureMap

Calculate the spatial probability density

\[\text{P}(r) = |\Psi_n(r)|^2\]

for each position \(r\) in system.positions where \(\Psi_n(r)\) is eigenvectors[:, n].

Parameters:
nint or array_like

Index of the desired eigenstate. If an array of indices is given, the probability will be calculated at each one and a sum will be returned.

reducefloat, optional

Reduce degenerate states by summing their probabilities. Neighboring states are considered degenerate if their energy is difference is lower than the value of reduce. This is disabled by passing reduce=0.

Returns:
StructureMap
calc_spatial_ldos(energy: float, broadening: float) StructureMap

Calculate the spatial local density of states at the given energy

\[\text{LDOS}(r) = \frac{1}{c \sqrt{2\pi}} \sum_n{|\Psi_n(r)|^2 e^{-\frac{(E_n - E)^2}{2 c^2}}}\]

for each position \(r\) in system.positions, where \(E\) is energy, \(c\) is broadening, \(E_n\) is eigenvalues[n] and \(\Psi_n(r)\) is eigenvectors[:, n].

Parameters:
energyfloat

The energy value for which the spatial LDOS is calculated.

broadeningfloat

Controls the width of the Gaussian broadening applied to the DOS.

Returns:
StructureMap
calc_wavefunction(k0: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], k1: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], *ks: Iterable[_SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes]], step: float = 0.1, point_labels: List[str] | None = None) Wavefunction

Calculate the wavefunction on a path in reciprocal space

Parameters:
k0, k1, *ksarray_like

Points in reciprocal space which form the path for the band calculation. At least two points are required.

stepfloat, optional

Calculation step length in reciprocal space units. Lower step values will return more detailed results.

point_labelsList[str], optional

The point_labels for plots

Returns:
Wavefunction
calc_wavefunction_area(k_area: Area) WavefunctionArea

Calculate the wavefunction on a path in reciprocal space

Parameters:
k_path~pybinding.Path

Points in reciprocal space which form the path for the band calculation. At least two points are required.

Returns:
WavefunctionArea
calc_wavefunction_path(k_path: Path) Wavefunction

Calculate the wavefunction on a path in reciprocal space

Parameters:
k_path~pybinding.Path

Points in reciprocal space which form the path for the band calculation. At least two points are required.

Returns:
Wavefunction
clear()

Clear the computed results and start over

static find_degenerate_states(energies: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], abs_tolerance: float = 1e-05) List[List[float]]

Return groups of indices which belong to degenerate states

Parameters:
energiesarray_like
abs_tolerancefloat, optional

Examples

>>> energies = np.array([0.1, 0.1, 0.2, 0.5, 0.5, 0.5, 0.7, 0.8, 0.8])
>>> Solver.find_degenerate_states(energies)
[[0, 1], [3, 4, 5], [7, 8]]
>>> energies = np.array([0.1, 0.2, 0.5, 0.7])
>>> Solver.find_degenerate_states(energies)
[]
report(shortform: bool = False) str

Return a report of the last solve() computation

Parameters:
shortformbool, optional

Return a short one line version of the report

set_wave_vector(k: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes])

Set the wave vector for periodic models

Parameters:
karray_like

Wave vector in reciprocal space.

solve()

Explicitly solve the eigenvalue problem right now

This method is usually not needed because the main result properties, eigenvalues and eigenvectors, will call this implicitly the first time they are accessed. However, since the solve() routine may be computationally expensive, it is useful to have the ability to call it ahead of time as needed.

property eigenvalues: ndarray

1D array of computed energy states

property eigenvectors: ndarray

2D array where each column represents a wave function

eigenvectors.shape == (system.num_sites, eigenvalues.size)

property get_wave_vector

Get the wave vector for periodic models

property model: Model

The tight-binding model attached to this solver

property system: System

The tight-binding system attached to this solver (shortcut for Solver.model.system)

arpack(model: Model, k: int, sigma: float = 0, **kwargs) Solver

ARPACK Solver implementation for sparse matrices

This solver is intended for large models with sparse Hamiltonian matrices. It only computes a small targeted subset of eigenvalues and eigenvectors. Internally this solver uses the scipy.sparse.linalg.eigsh() function for sparse Hermitian matrices.

Parameters:
modelModel

Model which will provide the Hamiltonian matrix.

kint

The desired number of eigenvalues and eigenvectors. This number must be smaller than the size of the matrix, preferably much smaller for optimal performance. The computed eigenvalues are the ones closest to sigma.

sigmafloat, optional

Look for eigenvalues near sigma.

**kwargs

Advanced arguments: forwarded to scipy.sparse.linalg.eigsh().

Returns:
Solver
dacp(model: Model, window: Tuple[float, float] = (-2, 2), random_vectors: int = 100, filter_order: int = 30, tol: float = 0.001, **kwargs) Solver | None

pyDACP Solver implementation for DACP method matrices

Some more text about DACP blablabla… Also look at https://gitlab.kwant-project.org/qt/pyDACP and install it with >>pip install git+https://gitlab.kwant-project.org/qt/pyDACP

Parameters:
modelpb.Model

Model which will provide the Hamiltonian matrix.

windowtuple of float

The lowest and highest eigenvalue between which to compute the solutions.

random_vectorsint, optional

Number of random vectors to use for the initial guess.

filter_orderint, optional

Order of the filter to use for the solution.

tolfloat, optional

Tolerance for the solution.

**kwargs

Advanced arguments: forwarded to scipy.sparse.linalg.eigsh().

Returns:
Solver
feast(model: Model, energy_range: Tuple[float, float], initial_size_guess: int, recycle_subspace: bool = False, is_verbose: bool = False) Solver

FEAST Solver implementation for sparse matrices

This solver is only available if the C++ extension module was compiled with FEAST.

Parameters:
modelModel

Model which will provide the Hamiltonian matrix.

energy_rangetuple of float

The lowest and highest eigenvalue between which to compute the solutions.

initial_size_guessint

Initial user guess for number of eigenvalues which will be found in the given energy_range. This value may be completely wrong - the solver will auto-correct as needed. However, for optimal performance the estimate should be as close to 1.5 * actual_size as possible.

recycle_subspacebool, optional

Reuse previously computed values as a starting point for the next computation. This improves performance when subsequent computations differ only slightly, as is the case for the band structure of periodic systems where the results change gradually as a function of the wave vector. It may hurt performance otherwise.

is_verbosebool, optional

Show the raw output from the FEAST routine.

Returns:
Solver
lapack(model: Model, **kwargs) Solver

LAPACK Solver implementation for dense matrices

This solver is intended for small models which are best represented by dense matrices. Always solves for all the eigenvalues and eigenvectors. Internally this solver uses the scipy.linalg.eigh() function for dense Hermitian matrices.

Parameters:
modelModel

Model which will provide the Hamiltonian matrix.

**kwargs

Advanced arguments: forwarded to scipy.linalg.eigh().

Returns:
Solver