sisl_toolbox.btd.DeviceGreen
- class sisl_toolbox.btd.DeviceGreen
Bases:
objectBlock-tri-diagonal Green function calculator
This class enables the extraction and calculation of some important quantities not currently accessible in TBtrans.
For instance it may be used to calculate scattering states from the Green function. Once scattering states have been calculated one may also calculate the eigenchannels.
Both calculations are very efficient and uses very little memory compared to the full matrices normally used.
Consider a regular 2 electrode setup with transport direction along the 3rd lattice vector. Then the following example may be used to calculate the eigen-channels.
The below short-form of reading all variables should cover most variables encountered in the FDF file.
G = DeviceGreen.from_fdf("RUN.fdf") # Calculate the scattering state from the left electrode # and then the eigen channels to the right electrode state = G.scattering_state("Left", E=0.1) eig_channel = G.eigenchannel(state, "Right")
The above
DeviceGreen.from_fdfis a short-hand for something like the below (it actually does more than that, so prefer thefrom_fdf):import sisl from sisl_toolbox.btd import * # First read in the required data H_elec = sisl.Hamiltonian.read("ELECTRODE.nc") H = sisl.Hamiltonian.read("DEVICE.nc") # remove couplings along the self-energy direction # to ensure no fake couplings. H.set_nsc(c=1) # Read in a single tbtrans output which contains the BTD matrices # and instructs this class how it should pivot the matrix to obtain # a BTD matrix. tbt = sisl.get_sile("siesta.TBT.nc") # Define the self-energy calculators which will downfold the # self-energies into the device region. # Since a downfolding will be done it requires the device Hamiltonian. H_elec.shift(tbt.mu("Left")) left = DownfoldSelfEnergy("Left", s.RecursiveSI(H_elec, "-C", eta=tbt.eta("Left"), tbt, H) H_elec.shift(tbt.mu("Right")) left = DownfoldSelfEnergy("Right", s.RecursiveSI(H_elec, "+C", eta=tbt.eta("Right"), tbt, H) G = DeviceGreen(H, [left, right], tbt)
Notes
Sometimes one wishes to investigate more details in the calculation process to discern importance of the eigenvalue separations.
When calculating scattering states/matrices one can reduce the complexity by removing eigen/singular values. By default we use the cutoff values as a relative cutoff value for the values. I.e. keeping
value / value.max() > cutoff. However, sometimes the relative value is a bad metric since there are still important values close to unity value. Consider e.g. an array of values of[1e5, 1e4, 1, 0.5, 1e-4]. In this case we would require the[1, 0.5]values as important, but this would only be grabbed by a relative cutoff value of1e-6which in some other cases are a too high value.Instead of providing cutoff values as float values, one can also pass a function that takes in an array of values. It should return the indices of the values it wishes to retain.
The below is equivalent to a cutoff value of
1e-4, or values above 0.01. >>> def cutoff_func(V): >>> return np.logical_or(V / V.max() > 1e-4, V > 1e-2).nonzero()[0]Passing functions for cutting off values can be useful because one can also debug the values and see what’s happening.
Methods
Hk(*args, **kwargs)Return the Hamiltonian matrix in the pivoted device region
Sk(*args, **kwargs)Return the overlap matrix in the pivoted device region
clear(*keys)Clean any memory used by this object
eigenchannel(state[, elec_to, ret_coeff])Calculate the eigenchannel from scattering states entering electrodes elec_to
eigenchannel_from_scattering_matrix(scat_matrix)Calculate the eigenchannel from a scattering matrix
from_fdf(fdf[, prefix, use_tbt_se, eta])Return a new
DeviceGreenusing information gathered from the fdf file.green(E[, k, format, dtype])Calculate the Green function for a given E and k point
scattering_matrix(E, elec_from[, elec_to, ...])Calculate the scattering matrix (S-matrix) between elec_from and elec_to
scattering_state(E, elec[, k, cutoff, ...])Calculate the scattering states for a given E and k point from a given electrode
spectral(E, elec[, k, format, method, herm, ...])Calculate the spectral function for a given E and k point from a given electrode
transmission(E, elec_from[, elec_to, k, dtype])Calculate the transmission between an electrode, and one or more other electrodes
- eigenchannel(state, elec_to=None, ret_coeff=False)[source]
Calculate the eigenchannel from scattering states entering electrodes elec_to
The energy and k-point is inferred from the state object as returned from
scattering_state.The eigenchannels are the eigenstates of the transmission matrix in the DOS weighted scattering states:
\[\begin{split}\mathbf A_{\mathfrak e_{\mathrm{from}} }(E,\mathbf k) \mathbf u_i &= 2\pi a_i \mathbf u_i \\ \mathbf t_{\mathbf u} &= \sum \langle \mathbf u | \boldsymbol\Gamma_{ \mathfrak e_{\mathrm{to}} } | \mathbf u\rangle\end{split}\]where the eigenvectors of \(\mathbf t_{\mathbf u}\) are the coefficients of the DOS weighted scattering states (\(\sqrt{2\pi a_i} u_i\)) for the individual eigen channels. The eigenvalues are the transmission eigenvalues. Further details may be found in [8].
- Parameters:
state (StateCElectron) – the scattering states as obtained from
scattering_stateelec_to (
strorint (listornot)) – which electrodes to consider for the transmission eigenchannel decomposition (the sum in the above equation). Defaults to all but the origin electrode.ret_coeff (bool) – return also the scattering state coefficients
- Returns:
T_eig (
sisl.physics.electron.StateCElectron) – the transmission eigenchannels, theT_eig.ccontains the transmission eigenvalues.coeff (
sisl.physics.electron.StateElectron) – coefficients of state that creates the transmission eigenchannels Only returned if ret_coeff is True. There is a one-to-one correspondance betweencoeffandT_eig(with a prefactor of \(\sqrt{2\pi}\)). This is equivalent to theT_eigstates in the scattering state basis.
- Return type:
- eigenchannel_from_scattering_matrix(scat_matrix, ret_out=False)[source]
Calculate the eigenchannel from a scattering matrix
The energy and k-point is inferred from the state_matrix object as returned from
scattering_matrix.The eigenchannels are the SVD of the scattering matrix in the DOS weighted scattering states:
- Parameters:
state_matrix – the scattering matrix as obtained from
scattering_matrixscat_matrix (StateElectron)
ret_out (bool)
- Returns:
T_eig_in (
sisl.physics.electron.StateCElectron) – the transmission eigenchannels as seen from the incoming state, theT_eig.ccontains the transmission eigenvalues.T_eig_out (
sisl.physics.electron.StateCElectron) – the transmission eigenchannels as seen from the outgoing state, theT_eig.ccontains the transmission eigenvalues. Only returned if ret_out is true.
- Return type:
- classmethod from_fdf(fdf, prefix='TBT', use_tbt_se=False, eta=None, **kwargs)[source]
Return a new
DeviceGreenusing information gathered from the fdf file.- Parameters:
fdf (SileLike) – fdf file to read the parameters from
prefix (Literal['TBT', 'TS']) – which prefix to use, if TBT it will prefer TBT prefix, but fall back to TS prefixes. If TS, only these prefixes will be used.
use_tbt_se (bool) – whether to use the TBT.SE.nc files for self-energies or calculate them on the fly.
eta (Optional[float]) – force a specific eta value
kwargs – passed to the class instantiating.
- Return type:
Self
- green(E, k=(0, 0, 0), format='array', dtype=np.complex128, **kwargs)[source]
Calculate the Green function for a given E and k point
The Green function is calculated as:
\[\mathbf G(E,\mathbf k) = \big[\mathbf S(\mathbf k) E - \mathbf H(\mathbf k) - \sum \boldsymbol \Sigma(E,\mathbf k)\big]^{-1}\]- Parameters:
E (complex) – the energy to calculate at, may be a complex value.
k (Sequence[float]) – k-point to calculate the Green function at
format (
{"array", "btd", "bm", "bd", "sparse"}) –return the matrix in a specific format
array: a regular numpy array (full matrix)
btd: a block-matrix object with only the diagonals and first off-diagonals
bm: a block-matrix object with diagonals and all off-diagonals
bd: a block-matrix object with only diagonals (no off-diagonals)
sparse: a sparse-csr matrix for the sparse elements as found in the Hamiltonian
dtype – the data-type of the array.
- Returns:
the Green function matrix, the format depends on format.
- Return type:
- scattering_matrix(E, elec_from, elec_to=None, k=(0, 0, 0), cutoff=1e-4, dtype=np.complex128)[source]
Calculate the scattering matrix (S-matrix) between elec_from and elec_to
The scattering matrix is calculated as
\[\mathbf S_{\mathfrak e'\mathfrak e}(E, \mathbf) = -\delta_{\alpha\beta} + i \tilde{\boldsymbol\Gamma}_{\mathfrak e'} \mathbf G \tilde{\boldsymbol\Gamma}_{\mathfrak e}\]Here \(\tilde{\boldsymbol\Gamma}\) is defined as:
\[\begin{split}\boldsymbol\Gamma(E,\mathbf k) \mathbf u_i &= \lambda_i \mathbf u_i \\ \tilde{\boldsymbol\Gamma}(E,\mathbf k) &= \operatorname{diag}\{ \sqrt{\boldsymbol\lambda} \} \mathbf u\end{split}\]Once the scattering matrices have been calculated one can calculate the full transmission function
\[\mathcal T_{\mathfrak e\to\mathfrak e'}(E, \mathbf k) = \operatorname{Tr}\big[ \mathbf S_{\mathfrak e'\mathfrak e }^\dagger \mathbf S_{\mathfrak e'\mathfrak e }\big]\]The scattering matrix approach can be found in details in [12].
- Parameters:
E (complex) – the energy to calculate at, may be a complex value.
elec_from (
strorint) – the electrode where the scattering matrix originates fromelec_to (
strorintorlist of) – where the scattering matrix ends in.k (Sequence[float]) – k-point to calculate the scattering matrix at
cutoff (float) – cutoff eigen states of the broadening matrix. The cutoff is based on a relative fraction of the maximum eigen value. that are below this value. For example, we keep according to \(\lambda_i/\max(\lambda_i) > \mathrm{cutoff}\). A too high value will remove too many eigen states and results will be wrong. A small value improves precision at the cost of bigger matrices. The \(\Gamma\) matrix should be positive definite, however, due to the imaginary part of the self-energies it tends to only be close to positive definite.
- Returns:
for each elec_to a scattering matrix will be returned. Its dimensions will be depending on the cutoff value at the cost of precision.
- Return type:
sisl.physics.electron.StateElectronortuple[sisl.physics.electron.StateElectron,...]
- scattering_state(E, elec, k=(0, 0, 0), cutoff=0.0, method='svd:gamma', dtype=np.complex128, *args, **kwargs)[source]
Calculate the scattering states for a given E and k point from a given electrode
The scattering states are the eigen states of the spectral function:
\[\mathbf A_{\mathfrak e}(E,\mathbf k) \mathbf u_i = 2\pi a_i \mathbf u_i\]where \(a_i\) is the DOS carried by the \(i\)’th scattering state.
- Parameters:
E (complex) – the energy to calculate at, may be a complex value.
elec (
strorint) – the electrode to calculate the spectral function fromk (Sequence[float]) – k-point to calculate the spectral function at
cutoff (float | Callable) – Cut off the returned scattering states at some DOS value. Any scattering states with relative eigenvalues (to the largest eigenvalue), lower than cutoff are discarded. For example, we keep according to \(\epsilon_i/\max(\epsilon_i) > \mathrm{cutoff}\). Values above or close to 0.1 should be used with care. Can be a function, see the details of this class.
method (Literal['svd:gamma', 'svd:A', 'eig']) – which method to use for calculating the scattering states. Use only the
eigmethod for testing purposes as it is extremely slow and requires a substantial amount of memory. Thesvd:gammais the fastests while retaining complete precision. Thesvd:Amay be even faster for very large systems with very little loss of precision, since it diagonalizes \(\mathbf A\) in the subspace of the electrode elec and reduces the propagated part of the spectral matrix.cutoff_elec (
float, optional) – Only used formethod=svd:A. The initial block of the spectral function is diagonalized and only eigenvectors with relative eigenvalues>=cutoff_elecare retained. thus reducing the initial propagated modes. The normalization explained for cutoff also applies here. Can be a function, see the details of this class.
- Returns:
the scattering states from the spectral function. The
scat.statecontains the scattering state vectors (eigenvectors of the spectral function).scat.ccontains the DOS of the scattering states scaled by \(1/(2\pi)\) so ensure correct density of states. One may recreate the spectral function withscat.outer(matrix=scat.c * 2 * pi).- Return type:
- spectral(E, elec, k=(0, 0, 0), format='array', method='column', herm=True, dtype=np.complex128, **kwargs)[source]
Calculate the spectral function for a given E and k point from a given electrode
The spectral function is calculated as:
\[\mathbf A_{\mathfrak{e}}(E,\mathbf k) = \mathbf G(E,\mathbf k)\boldsymbol\Gamma_{\mathfrak{e}}(E,\mathbf k) \mathbf G^\dagger(E,\mathbf k)\]- Parameters:
E (complex) – the energy to calculate at, may be a complex value.
elec (
strorint) – the electrode to calculate the spectral function fromk (Sequence[float]) – k-point to calculate the spectral function at
format (
{"array", "btd", "bm", "bd"}) –return the matrix in a specific format
array: a regular numpy array (full matrix)
bm: in block-matrix form (full matrix)
btd: a block-matrix object with only the diagonals and first off-diagonals
bd: same as btd, since the off-diagonals are already calculated
method (Literal['column', 'propagate']) – which method to use for calculating the spectral function. Depending on the size of the BTD blocks one may be faster than the other. For large systems you are recommended to time the different methods and stick with the fastest one, they are numerically identical.
herm (bool) – The spectral function is a Hermitian matrix, by default (True), the methods that can utilize the Hermitian property only calculates the lower triangular part of \(\mathbf A\), and then copies the Hermitian to the upper part. By setting this to False the entire matrix is explicitly calculated.
- Returns:
the spectral function for a given electrode in the format as specified by format. Note that some formats does not calculate the entire spectral function matrix.
- Return type:
- transmission(E, elec_from, elec_to=None, k=(0, 0, 0), dtype=np.complex128)[source]
Calculate the transmission between an electrode, and one or more other electrodes
The transmission function is calculated as:
\[\mathcal T_{\mathfrak e\to\mathfrak e'}(E,\mathbf k) = \boldsymbol \Gamma_{\mathfrak e'}(E,\mathbf k) \mathbf G(E,\mathbf k) \boldsymbol\Gamma_{\mathfrak e}(E,\mathbf k) \mathbf G^\dagger(E,\mathbf k)\]- Parameters:
- Return type: