3. kwantspectrum
– Reference documentation¶
-
kwantspectrum.
spectrum
(syst, args=(), *, params=None, kmin=-3.141592653589793, kmax=3.141592653589793, orderpoint=0, tol=1e-08, match=<function _match_functions>)¶ Interpolate the dispersion function and provide methods to simplify curve sketching and analyzation the periodic spectrum.
- Parameters
syst (kwant.system.InfiniteSystem) – The low level infinite system for which the energies are to be calculated.
args (tuple, defaults to empty) – Positional arguments to pass to the
hamiltonian
method. Mutually exclusive with ‘params’.params (dict, optional) – Dictionary of parameter names and their values. Mutually exclusive with ‘args’.
kmin (float, optional) – Left-hand site of the momentum interval over which the band structure should be analyzed. Defaut = \(-\pi\)
kmax (float, optional) – Right-hand site of the momentum interval over which the band structure should be analyzed. Defaut = \(\pi\)
orderpoint (int or float, optional) – Momentum value where the band ordering is defined. Bands are ordered in energy, such that the band with the lowest energy at momentum orderpoint has index 0. Defaut = 0
tol (float, optional) – Numerical tolerance of the interpolated function.
match (callable, optional) – Matching algorithm
- Returns
spec
- Return type
BandSketching
instance or a list of these
Notes
The tolerance \(tol\) is the required precision of the interpolated band structure. The matching algorithm which orders the bands should work correctly quite independently of the tolerance \(tol\). However, as the matching algorithm uses piecewise cubic interpolations for the error estimate internally, the cubic piecewise interpolation should be used for the subsequent interpolation of the spectrum for consistency.
-
kwantspectrum.
spectra
(syst, *args, **kwargs)¶ Interpolates a sequence of dispersion functions and provide methods to simplify curve sketching and analyzation the periodic spectra.
- Parameters
syst (
kwant.system.InfiniteSystem
or sequence thereof) – The low level infinite systems.- Returns
specs
- Return type
list of
BandSketching
Notes
The function is similar to
spectrum
but works as well for a sequence of leads. It always returns a list of spectra. Seespectrum
for optional function arguments.
-
kwantspectrum.
intersect_intervals
(interval_a, interval_b, tol=1e-16)¶ Return the intersecting intervals between two lists of intervals.
-
class
kwantspectrum.
BandSketching
(x, y, dy, mode_function, tol=1e-08, interpolation=<function _cubic_interpolation>)¶ Interpolate the dispersion function and provide methods to simplify curve sketching and analyzation the periodic spectrum.
- Parameters
x (numpy list, shape (n, )) – momentum (k) points
y (numpy list, shape (n, nbands)) – energy values
dy (numpy list, shape (n, nbands)) – velocity values
mode_function (callable) – calculates the open modes for a given energy
tol (float) – Numerical accuracy of the interpolated function.
interpolation (object, optional) – Method to perform numerical interpolations. Default = cubic Hermite interpolation (piecewise cubic interpolation, 0. and 1. order of the exact and interpolation function are similar on the knots).
Notes
The routine designed for a periodic function that is represented as discrete values on a finite gridpoints. Finite accuracy of the data and possible unphysical results are taking into account by rounding value \(tol\).
-
energy_to_scattering_mode
(energy, band, kmin, kmax)¶ Finds the scattering mode index for a band giving its energy.
- Parameters
energy (float) –
band (int) – band index
kmax (kmin,) – momentum interval
[kmin, kmax]
including the searched momentum, where the velocity does not change sign
- Returns
mode – kwant scattering mode index The mode index fulfills: 0 <= mode < nbands where nbands is the total number of bands of the spectrum. If no open mode could be found, mode = -1 is returned.
- Return type
int
Notes
An exception is raised if the momentum interval [kmin, kmax] is badly chose, such that either no mode / multiple modes are found, and therefore no unique (energy, band) -> mode mapping is possible, a warning is printed.
-
intersect
(f, band, derivative_order=0, kmin=None, kmax=None, tol=None, ytol=None)¶ Returns all momentum (k) points, that solves the equation: \(\partial_k^{n} E(k) = f(k),\, k_{min} \leq k \leq k_{max}\).
- Parameters
f (scalar numerical value or callable) – Equation to solve \(\partial_k^{n} E(k) = f\). If f is callable, solve equation: \(\partial_k^{n} E(k) = f(k)\).
band (int) – band index, requirement: 0 <= band < nbands.
derivative_order (int, optional) – Derivative order (n) of the band dispersion. Default is zero.
kmin (scalar numeric value, optional) – Lowest k point value. Default is kmin from initialization.
kmax (scalar numeric value, optional) – Largest k point value. Default is kmax from initialization.
tol (float, optional) – Numerical tolerance, k points closer tol are merged to the same point. Default is the tol from initialization.
ytol (float, optional) – Numerical tolerance to remove noise if the spectrum \(\partial_k^{n} E(k)\) is almost flat. Values for the spectrum are set to thier mean value (averaged over all momentum points where the band is sampled), if they flucutate more than ytol. Default is the tol from initialization.
- Returns
k – List of momentum (k) points, that solves the above equation.
- Return type
numpy list.
Notes
kmin and kmin can be larger than the first Brillouin zone. In that case we return also the periodic images.
-
intervals
(band, derivative_order=0, lower=None, upper=None, kmin=None, kmax=None, tol=None)¶ returns a list of momentum (k) intervals that solves the equation \(\text{lower} \leq \partial^m_k E_n(k) \leq \text{upper}\)
- Parameters
band (int) – band index n, requirement: 0 <= n < nbands.
derivative_order (int, optional) – Derivative order m of the band dispersion. Default is zero.
lower (int or float, optional) – Select intervals such that \(\partial_k^{m} E_n(k) \geq \text{lower}\)
upper (int or float, optional) – Select intervals such that \(\partial_k^{m} E_n(k) \leq \text{upper}\)
kmin (scalar numeric value, optional) – Lowest k point value. Default is kmin from initialization.
kmax (scalar numeric value, optional) – Largest k point value. Default is kmax from initialization.
tol (float, optional) – Tolerance of the interval selection, such that intervals smaller tol are neglected and two intervals closer tol are merged.
- Returns
intervals – Ordered list of momentum intervals \([k_i, k_{i+1}]\). Each interval is given in form of a tuple.
- Return type
list.
-
momentum_to_scattering_mode
(k, band)¶ Finds the scattering mode index for a band with momentum k.
- Parameters
k (float) – momentum value
band (int) – band index
- Returns
mode – kwant scattering mode index. The mode index fulfills: 0 <= mode < nbands where nbands is the total number of bands of the spectrum. If no open mode could be found, mode = -1 is returned.
- Return type
int
-
set_period
(period)¶ Set periodicity of the spectrum
- Parameters
period (2-tuple or False) – Perodic interval provided in the form (left_bound, right_bound) over which the band structure is assumed as periodic. If period is set to false, the dispersion function is evaluated as non-periodic fuction and the user has to perform periodic mapping herself if needed.