chebyshev

Computations based on Chebyshev polynomial expansion

The kernel polynomial method (KPM) can be used to approximate various functions by expanding them in a series of Chebyshev polynomials.

Classes

KPM(impl)

The common interface for various KPM implementations

_ComputeProgressReporter()

_PythonImpl(model, energy_range, kernel, **_)

Basic Python/SciPy implementation of KPM

Functions

_kpm_python(model[, energy_range, kernel])

Basic Python/SciPy implementation of KPM

dirichlet_kernel()

The Dirichlet kernel -- returns raw moments, least favorable choice

jackson_kernel()

The Jackson kernel -- a good general-purpose kernel, appropriate for most applications

kpm(model[, energy_range, kernel, ...])

The default CPU implementation of the Kernel Polynomial Method

kpm_cuda(model[, energy_range, kernel])

Same as kpm() except that it's executed on the GPU using CUDA (if supported)

lorentz_kernel([lambda_value])

The Lorentz kernel -- best for Green's function

class KPM(impl: kpm)

The common interface for various KPM implementations

It should not be created directly but via specific functions like kpm() or kpm_cuda().

All implementations are based on: https://doi.org/10.1103/RevModPhys.78.275

__call__(*args, **kwargs) ndarray

Call self as a function.

calc_conductivity(chemical_potential: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], broadening: float, temperature: float, direction: Literal['xx', 'xy', 'xz', 'yx', 'yy', 'yz', 'zx', 'zy', 'zz'] = 'xx', volume: float = 1.0, num_random: int = 1, num_points: int = 1000) Series

Calculate Kubo-Bastin electrical conductivity as a function of chemical potential

The return value is in units of the conductance quantum (e^2 / hbar) not taking into account spin or any other degeneracy.

The calculation is based on: https://doi.org/10.1103/PhysRevLett.114.116602.

Parameters:
chemical_potentialarray_like

Values (in eV) for which the conductivity is calculated.

broadeningfloat

Width (in eV) of the smallest detail which can be resolved in the chemical potential. Lower values result in longer calculation time.

temperaturefloat

Value of temperature for the Fermi-Dirac distribution.

directionOptional[str]

Direction in which the conductivity is calculated. E.g., “xx”, “xy”, “zz”, etc.

volumeOptional[float]

The volume of the system.

num_randomint

The number of random vectors to use for the stochastic calculation of KPM moments. Larger numbers improve the quality of the result but also increase calculation time linearly. Fortunately, result quality also improves with system size, so the DOS of very large systems can be calculated accurately with only a small number of random vectors.

num_pointsOptional[int]

Number of points for integration.

Returns:
Series
calc_dos(energy: ndarray, broadening: float, num_random: int = 1) Series

Calculate the density of states as a function of energy

Parameters:
energyndarray

Values for which the DOS is calculated.

broadeningfloat

Width, in energy, of the smallest detail which can be resolved. Lower values result in longer calculation time.

num_randomint

The number of random vectors to use for the stochastic calculation of KPM moments. Larger numbers improve the quality of the result but also increase calculation time linearly. Fortunately, result quality also improves with system size, so the DOS of very large systems can be calculated accurately with only a small number of random vectors.

Returns:
Series
calc_greens(i: int, j: int, energy: ndarray, broadening: float) ndarray

Calculate Green’s function of a single Hamiltonian element

Parameters:
iint or list

Hamiltonian index.

jint

Hamiltonian index.

energyndarray

Energy value array.

broadeningfloat

Width, in energy, of the smallest detail which can be resolved. Lower values result in longer calculation time.

Returns:
ndarray

Array of the same size as the input energy.

calc_ldos(energy: ndarray, 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

Parameters:
energyndarray

Values for which the LDOS is calculated.

broadeningfloat

Width, in energy, of the smallest detail which can be resolved. Lower values result in longer calculation time.

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 be 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_spatial_ldos(energy: ndarray, broadening: float, shape: Shape, sublattice: str = '') SpatialLDOS

Calculate the LDOS as a function of energy and space (in the area of the given shape)

Parameters:
energyndarray

Values for which the LDOS is calculated.

broadeningfloat

Width, in energy, of the smallest detail which can be resolved. Lower values result in longer calculation time.

shapeShape

Determines the site positions at which to do the calculation.

sublatticestr

Only look for sites of a specific sublattice, within the shape. The default value considers any sublattice.

Returns:
SpatialLDOS
deferred_ldos(energy: ndarray, 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 = '') DeferredXd

Same as calc_ldos() but for parallel computation: see the parallel module

Parameters:
energyndarray

Values for which the LDOS is calculated.

broadeningfloat

Width, in energy, of the smallest detail which can be resolved. Lower values result in longer calculation time.

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.

Returns:
moments(num_moments: int, alpha: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], beta: _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes] | None = None, op: csr_matrix | None = None) ndarray

Calculate KPM moments in the form of expectation values

The result is an array of moments where each value is equal to:

\[\mu_n = <\beta|op \cdot T_n(H)|\alpha>\]
Parameters:
num_momentsint

The number of moments to calculate.

alphaarray_like

The starting state vector of the KPM iteration.

betaOptional[array_like]

If not given, defaults to \(\beta = \alpha\).

opOptional[csr_matrix]

Operator in the form of a sparse matrix. If omitted, an identity matrix is assumed: \(\mu_n = <\beta|T_n(H)|\alpha>\).

Returns:
ndarray
report(shortform: bool = False) str

Return a report of the last computation

Parameters:
shortformbool, optional

Return a short one line version of the report

property block_diagonal: List[int]

The first index of the reordered matrix where a block of a block-diagonal matrix ends.

property kernel: KPMKernel

The damping kernel

property model: Model

The tight-binding model holding the Hamiltonian

property scaling_factors: Tuple[float, float]

A tuple of KPM scaling factors a and b

property system: System

The tight-binding system (shortcut for KPM.model.system)

property zero_row: List[int]

The index of a row of zeros in the reordered matrix.

dirichlet_kernel() KPMKernel

The Dirichlet kernel – returns raw moments, least favorable choice

This kernel doesn’t modify the moments at all. The resulting moments represent just a truncated series which results in lots of oscillation in the reconstructed function. Therefore, this kernel should almost never be used. It’s only here in case the raw moment values are needed for some other purpose. Note that required_num_moments() returns N = pi / sigma for compatibility with the Jackson kernel, but there is no actual broadening associated with the Dirichlet kernel.

jackson_kernel() KPMKernel

The Jackson kernel – a good general-purpose kernel, appropriate for most applications

Imposes Gaussian broadening sigma = pi / N where N is the number of moments. The broadening value is user-defined for each function calculation (LDOS, Green’s, etc.). The number of moments is then determined based on the broadening – it’s not directly set by the user.

kpm(model: Model, energy_range: Tuple[float, float] | None = None, kernel: KPMKernel | Literal['default'] = 'default', num_threads: int | Literal['auto'] = 'auto', silent: bool = False, **kwargs) KPM

The default CPU implementation of the Kernel Polynomial Method

This implementation works on any system and is well optimized.

Parameters:
modelModel

Model which will provide the Hamiltonian matrix.

energy_rangeOptional[Tuple[float, float]]

KPM needs to know the lowest and highest eigenvalue of the Hamiltonian, before computing the expansion moments. By default, this is determined automatically using a quick Lanczos procedure. To override the automatic boundaries pass a (min_value, max_value) tuple here. The values can be overestimated, but note that performance drops as the energy range becomes wider. On the other hand, underestimating the range will produce NaN values in the results.

kernelKernel

The kernel in the Kernel Polynomial Method. Used to improve the quality of the function reconstructed from the Chebyshev series. Possible values are jackson_kernel() or lorentz_kernel(). The Jackson kernel is used by default.

num_threadsint

The number of CPU threads to use for calculations. This is automatically set to the number of logical cores available on the current machine.

silentbool

Don’t show any progress messages.

Returns:
KPM
kpm_cuda(model: Model, energy_range: Tuple[float, float] | None = None, kernel: KPMKernel | Literal['default'] = 'default', **kwargs) KPM

Same as kpm() except that it’s executed on the GPU using CUDA (if supported)

See kpm() for detailed parameter documentation. This method is only available if the C++ extension module was compiled with CUDA.

Parameters:
modelModel
energy_rangeOptional[Tuple[float, float]]
kernelKernel
Returns:
KPM
lorentz_kernel(lambda_value: float = 4.0) KPMKernel

The Lorentz kernel – best for Green’s function

This kernel is most appropriate for the expansion of the Green’s function because it most closely mimics the divergences near the true eigenvalues of the Hamiltonian. The Lorentzian broadening is given by epsilon = lambda / N where N is the number of moments.

Parameters:
lambda_valuefloat

May be used to fine-tune the smoothness of the convergence. Usual values are between 3 and 5. Lower values will speed up the calculation at the cost of accuracy. If in doubt, leave it at the default value of 4.