sisl.io.tbtrans.tbtncSileTBtrans

class sisl.io.tbtrans.tbtncSileTBtrans(filename, mode='r', lvl=0, access=1, *args, **kwargs)

Bases: sisl.io.tbtrans._devncSileTBtrans

TBtrans output file object

Implementation of the TBtrans output *.TBT.nc files which contains calculated quantities related to the NEGF code TBtrans.

Although the TBtrans code is in fortran and the resulting NetCDF file variables are in fortran indexing (1-based), everything is returned as Python indexing (0-based) when using Python scripts.

In the following equations we will use this notation:

  • \(\alpha\) and \(\beta\) are atomic indices

  • \(\nu\) and \(\mu\) are orbital indices

A word on DOS normalization:

All the device region DOS functions may request a normalization depending on a variety of functions. You are highly encouraged to read the documentation for the norm function and to consider the benefit of using the norm='atom' normalization to more easily compare various partitions of DOS.

Notes

The API for this class are largely equivalent to the arguments of the sdata command-line tool, with the execption that the command-line tool uses Fortran indexing numbers (1-based).

Methods

ADOS([elec, E, kavg, atoms, orbitals, sum, norm])

Spectral density of states (DOS) (1/eV).

Adensity_matrix(elec, E[, kavg, isc, ...])

Spectral function density matrix at energy E (1/eV)

BDOS([elec, E, kavg, sum, norm])

Bulk density of states (DOS) (1/eV).

DOS([E, kavg, atoms, orbitals, sum, norm])

Green function density of states (DOS) (1/eV).

Eindex(E)

Return the closest energy index corresponding to the energy E

a2p(atoms)

Return the pivoting orbital indices (0-based) for the atoms, possibly on an electrode

a_down(elec[, bulk])

Down-folding atomic indices for a given electrode

a_elec(elec)

Electrode atomic indices for the full geometry (sorted)

atom_ACOHP(elec, E[, kavg, isc, orbitals, uc])

Atomic COHP curve of the spectral function

atom_ACOOP(elec, E[, kavg, isc, orbitals, uc])

Atomic COOP curve of the spectral function

atom_COHP(E[, kavg, isc, orbitals, uc])

Atomic COHP curve of the Green function

atom_COHP_from_orbital(COHP[, uc])

Calculate the atomic COHP curve from the orbital COHP

atom_COOP(E[, kavg, isc, orbitals, uc])

Atomic COOP curve of the Green function

atom_COOP_from_orbital(COOP[, uc])

Calculate the atomic COOP curve from the orbital COOP

atom_current(elec, E[, kavg, activity, orbitals])

Atomic current of atoms

atom_current_from_orbital(Jij[, activity])

Atomic current of atoms by passing the orbital current

bloch(elec)

Bloch-expansion coefficients for an electrode

bond_current(elec, E[, kavg, isc, only, ...])

Bond-current between atoms (sum of orbital currents)

bond_current_from_orbital(Jij[, only, uc])

Bond-current between atoms (sum of orbital currents) from an external orbital current

btd([elec])

Block-sizes for the BTD method in the device/electrode region

chemical_potential(elec)

Return the chemical potential associated with the electrode elec

close()

current([elec_from, elec_to, kavg])

Current from from to to using the k-weights and energy spacings in the file.

current_parameter(elec_from, mu_from, ...[, ...])

Current from from to to using the k-weights and energy spacings in the file.

density_matrix(E[, kavg, isc, orbitals, ...])

Density matrix from the Green function at energy E (1/eV)

dir_file([filename, filename_base])

File of the current Sile

electron_temperature(elec)

Electron bath temperature [Kelvin]

eta([elec])

The imaginary part used when calculating the self-energies in eV (or for the device

fano([elec_from, elec_to, kavg, zero_T])

The Fano-factor for the calculation (requires calculated transmission eigenvalues)

info([elec])

Information about the calculated quantities available for extracting in this file

iter([group, dimension, variable, levels, root])

Iterator on all groups, variables and dimensions.

kT(elec)

Electron bath temperature [eV]

kindex(k)

Return the index of the k-point that is closests to the queried k-point (in reduced coordinates)

mu(elec)

Return the chemical potential associated with the electrode elec

n_btd([elec])

Number of blocks in the BTD partioning

na_down(elec)

Number of atoms in the downfolding region (without device downfolded region)

no_down(elec)

Number of orbitals in the downfolding region (plus device downfolded region)

no_e(elec)

Number of orbitals in the downfolded region of the electrode in the device

noise_power([elec_from, elec_to, kavg])

Noise power from to to using the k-weights and energy spacings in the file (temperature dependent)

norm([atoms, orbitals, norm])

Normalization factor depending on the input

o2p(orbitals[, elec])

Return the pivoting indices (0-based) for the orbitals, possibly on an electrode

orbital_ACOHP(elec, E[, kavg, isc, orbitals])

Orbital resolved COHP analysis of the spectral function

orbital_ACOOP(elec, E[, kavg, isc, orbitals])

Orbital COOP analysis of the spectral function

orbital_COHP(E[, kavg, isc, orbitals])

Orbital resolved COHP analysis of the Green function

orbital_COOP(E[, kavg, isc, orbitals])

Orbital COOP analysis of the Green function

orbital_current(elec, E[, kavg, isc, only, ...])

Orbital current originating from elec as a sparse matrix

pivot([elec, in_device, sort])

Return the pivoting indices for a specific electrode (in the device region) or the device

pivot_down(elec)

Pivoting orbitals for the downfolding region of a given electrode

read(*args, **kwargs)

Generic read method which should be overloaded in child-classes

read_data(*args, **kwargs)

Read specific type of data.

read_geometry(*args, **kwargs)

Returns Geometry object from this file

read_supercell()

Returns SuperCell object from this file

reflection([elec, kavg, from_single])

Reflection into electrode elec

shot_noise([elec_from, elec_to, classical, kavg])

Shot-noise term from to to using the k-weights

transmission([elec_from, elec_to, kavg])

Transmission from elec_from to elec_to.

transmission_bulk([elec, kavg])

Bulk transmission for the elec electrode

transmission_eig([elec_from, elec_to, kavg])

Transmission eigenvalues from elec_from to elec_to.

vector_current(elec, E[, kavg, only, orbitals])

Vector for each atom describing the mean path for the current travelling through the atom

vector_current_from_bond(Jab)

Vector for each atom being the sum of bond-current times the normalized bond between the atoms

write(*args, **kwargs)

Generic write method which should be overloaded in child-classes

write_tbtav(*args, **kwargs)

Convert this to a TBT.AV.nc file, i.e. all k dependent quantites are averaged out.

E

Sampled energy-points in file

a_buf

Atomic indices (0-based) of device atoms

a_dev

Atomic indices (0-based) of device atoms (sorted)

base_file

File of the current Sile

cell

Unit cell in file

elecs

List of electrodes

file

File of the current Sile

geom

The associated geometry from this file

geometry

The associated geometry from this file

k

Sampled k-points in file

kpt

Sampled k-points in file

lasto

Last orbital of corresponding atom

nE

Number of energy-points in file

na

Returns number of atoms in the cell

na_b

Number of atoms in the buffer region

na_buffer

Number of atoms in the buffer region

na_d

Number of atoms in the device region

na_dev

Number of atoms in the device region

na_u

Returns number of atoms in the cell

ne

Number of energy-points in file

nk

Number of k-points in file

nkpt

Number of k-points in file

no

Returns number of orbitals in the cell

no_d

Number of orbitals in the device region

no_u

Returns number of orbitals in the cell

o_dev

Orbital indices (0-based) of device orbitals (sorted)

wk

Weights of k-points in file

wkpt

Weights of k-points in file

xa

Atomic coordinates in file

xyz

Atomic coordinates in file

ADOS(elec=0, E=None, kavg=True, atoms=None, orbitals=None, sum=True, norm='none')[source]

Spectral density of states (DOS) (1/eV).

Extract the spectral DOS from electrode elec on a selected subset of atoms/orbitals in the device region

\[\mathrm{ADOS}_\mathfrak{el}(E) = \frac{1}{2\pi N} \sum_{\nu\in \mathrm{atom}/\mathrm{orbital}} [\mathbf{G}(E)\Gamma_\mathfrak{el}\mathbf{G}^\dagger]_{\nu\nu}(E)\]

The normalization constant (\(N\)) is defined in the routine norm and depends on the arguments.

Parameters
  • elec (str, int, optional) – electrode originating spectral function

  • E (float or int, optional) – optionally only return the DOS of atoms at a given energy point

  • kavg (bool, int, optional) – whether the returned DOS is k-averaged, or an explicit (unweighed) k-point is returned

  • atoms (array_like of int or bool, optional) – only return for a given set of atoms (default to all). NOT allowed with orbitals keyword. If True it will use all atoms in the device. False is equivalent to None.

  • orbitals (array_like of int or bool, optional) – only return for a given set of orbitals (default to all) NOT allowed with atoms keyword. If True it will use all orbitals in the device. False is equivalent to None.

  • sum (bool, optional) – whether the returned quantities are summed or returned as is, i.e. resolved per atom/orbital.

  • norm ({'none', 'atom', 'orbital', 'all'}) – how the normalization of the summed DOS is performed (see norm routine).

See also

DOS

the total density of states (including bound states)

BDOS

the bulk density of states in an electrode

Adensity_matrix(elec, E, kavg=True, isc=None, orbitals=None, geometry=None)[source]

Spectral function density matrix at energy E (1/eV)

The density matrix can be used to calculate the LDOS in real-space.

The \(\mathrm{LDOS}(E, \mathbf r)\) may be calculated using the density routine. Basically the LDOS in real-space may be calculated as

\[\rho_{\mathbf A_{\mathfrak{el}}}(E, \mathbf r) = \frac{1}{2\pi}\sum_{\nu\mu}\phi_\nu(\mathbf r)\phi_\mu(\mathbf r) \Re[\mathbf A_{\mathfrak{el}, \nu\mu}(E)]\]

where \(\phi\) are the orbitals. Note that the broadening used in the TBtrans calculations ensures the broadening of the density, i.e. it should not be necessary to perform energy averages over the density matrices.

Parameters
  • elec (str or int) – the electrode of originating electrons

  • E (float or int) – the energy or the energy index of density matrix. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned density matrix is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned density matrix from unit-cell ([None, None, None]) to the given supercell, the default is all density matrix elements for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain density matrix elements for a subset of orbitals, all other are set to 0.

  • geometry (Geometry, optional) – geometry that will be associated with the density matrix. By default the geometry contained in this file will be used. However, then the atomic species are probably incorrect, nor will the orbitals contain the basis-set information required to generate the required density in real-space.

See also

density_matrix

Green function density matrix

Returns

object containing the Geometry and the density matrix elements

Return type

DensityMatrix

BDOS(elec=0, E=None, kavg=True, sum=True, norm='none')[source]

Bulk density of states (DOS) (1/eV).

Extract the bulk DOS from electrode elec on a selected subset of atoms/orbitals in the device region

\[\mathrm{BDOS}_\mathfrak{el}(E) = -\frac{1}{\pi} \Im\mathbf{G}(E)\]
Parameters
  • elec (str, int, optional) – electrode where the bulk DOS is returned

  • E (float or int, optional) – optionally only return the DOS of atoms at a given energy point

  • kavg (bool, int, optional) – whether the returned DOS is k-averaged, or an explicit (unweighed) k-point is returned

  • sum (bool, optional) – whether the returned quantities are summed or returned as is, i.e. resolved per atom/orbital.

  • norm ({'none', 'atom', 'orbital', 'all'}) – whether the returned quantities are summed or normed by total number of orbitals. Currently one cannot extract DOS per atom/orbital.

See also

DOS

the total density of states (including bound states)

ADOS

the spectral density of states from an electrode

DOS(E=None, kavg=True, atoms=None, orbitals=None, sum=True, norm='none')[source]

Green function density of states (DOS) (1/eV).

Extract the DOS on a selected subset of atoms/orbitals in the device region

\[\mathrm{DOS}(E) = -\frac{1}{\pi N} \sum_{\nu\in \mathrm{atom}/\mathrm{orbital}} \Im \mathbf{G}_{\nu\nu}(E)\]

The normalization constant (\(N\)) is defined in the routine norm and depends on the arguments.

Parameters
  • E (float or int, optional) – optionally only return the DOS of atoms at a given energy point

  • kavg (bool, int, optional) – whether the returned DOS is k-averaged, or an explicit (unweighed) k-point is returned

  • atoms (array_like of int or bool, optional) – only return for a given set of atoms (default to all). NOT allowed with orbitals keyword. If True it will use all atoms in the device. False is equivalent to None.

  • orbitals (array_like of int or bool, optional) – only return for a given set of orbitals (default to all) NOT allowed with atoms keyword. If True it will use all orbitals in the device. False is equivalent to None.

  • sum (bool, optional) – whether the returned quantities are summed or returned as is, i.e. resolved per atom/orbital.

  • norm ({'none', 'atom', 'orbital', 'all'}) – how the normalization of the summed DOS is performed (see norm routine)

See also

ADOS

the spectral density of states from an electrode

BDOS

the bulk density of states in an electrode

property E

Sampled energy-points in file

Eindex(E)

Return the closest energy index corresponding to the energy E

Parameters

E (float or int or str) – if int, return it-self, else return the energy index which is closests to the energy. For a str it will be parsed to a float and treated as such.

__init__(filename, mode='r', lvl=0, access=1, *args, **kwargs)
a2p(atoms)

Return the pivoting orbital indices (0-based) for the atoms, possibly on an electrode

This is equivalent to:

>>> p = self.o2p(self.geometry.a2o(atom, True))

Will warn if an atom requested is not in the device list of atoms.

Parameters

atoms (array_like or int) – atomic indices (0-based)

property a_buf

Atomic indices (0-based) of device atoms

property a_dev

Atomic indices (0-based) of device atoms (sorted)

a_down(elec, bulk=False)

Down-folding atomic indices for a given electrode

Parameters
  • elec (str or int) – electrode to retrieve indices for

  • bulk (bool, optional) – whether the returned indices are only in the pristine electrode, or the down-folding region (electrode + downfolding region, not in device)

a_elec(elec)

Electrode atomic indices for the full geometry (sorted)

Parameters

elec (str or int) – electrode to retrieve indices for

atom_ACOHP(elec, E, kavg=True, isc=None, orbitals=None, uc=False)[source]

Atomic COHP curve of the spectral function

Parameters
  • elec (str or int) – the electrode of the spectral function

  • E (float or int) – the energy or the energy index of COHP. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned COHP is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned COHP from unit-cell ([None, None, None]) to the given supercell, the default is all COHP for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain COHP matrix elements for a subset of orbitals, all other are set to 0.

  • uc (bool, optional) – whether the returned COHP are only in the unit-cell. If True this will return a sparse matrix of shape = (self.na, self.na), else, it will return a sparse matrix of shape = (self.na, self.na * self.n_s). One may figure out the connections via sc_index.

See also

orbital_COHP

orbital resolved COHP analysis of the Green function

atom_COHP_from_orbital

transfer an orbital COHP to atomic COHP

atom_COHP

atomic COHP analysis of the Green function

atom_ACOOP

atomic COOP analysis of the spectral function

atom_ACOOP(elec, E, kavg=True, isc=None, orbitals=None, uc=False)[source]

Atomic COOP curve of the spectral function

Parameters
  • elec (str or int) – the electrode of the spectral function

  • E (float or int) – the energy or the energy index of COOP. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned COOP is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned COOP from unit-cell ([None, None, None]) to the given supercell, the default is all COOP for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain COOP matrix elements for a subset of orbitals, all other are set to 0.

  • uc (bool, optional) – whether the returned COOP are only in the unit-cell. If True this will return a sparse matrix of shape = (self.na, self.na), else, it will return a sparse matrix of shape = (self.na, self.na * self.n_s). One may figure out the connections via sc_index.

See also

orbital_COOP

orbital resolved COOP analysis of the Green function

atom_COOP_from_orbital

transfer an orbital COOP to atomic COOP

atom_COOP

atomic COOP analysis of the Green function

atom_ACOHP

atomic COHP analysis of the spectral function

atom_COHP(E, kavg=True, isc=None, orbitals=None, uc=False)[source]

Atomic COHP curve of the Green function

Parameters
  • E (float or int) – the energy or the energy index of COHP. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned COHP is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned COHP from unit-cell ([None, None, None]) to the given supercell, the default is all COHP for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain COHP matrix elements for a subset of orbitals, all other are set to 0.

  • uc (bool, optional) – whether the returned COHP are only in the unit-cell. If True this will return a sparse matrix of shape = (self.na, self.na), else, it will return a sparse matrix of shape = (self.na, self.na * self.n_s). One may figure out the connections via sc_index.

See also

orbital_COHP

orbital resolved COHP analysis of the Green function

atom_COHP_from_orbital

transfer an orbital COHP to atomic COHP

atom_ACOHP

atomic COHP analysis of the spectral function

atom_COOP

atomic COOP analysis of the Green function

atom_COHP_from_orbital(COHP, uc=False)[source]

Calculate the atomic COHP curve from the orbital COHP

The atomic COHP are a sum over all orbital COHP:

\[\mathrm{COHP}_{\alpha\beta} = \sum_{\nu\in\alpha}\sum_{\mu\in\beta} \mathrm{COHP}_{\nu\mu}\]
Parameters
  • COHP (scipy.sparse.csr_matrix) – the orbital COHP as retrieved from orbital_COHP or orbital_ACOHP

  • uc (bool, optional) – whether the returned COHP are only in the unit-cell. If True this will return a sparse matrix of shape = (self.na, self.na), else, it will return a sparse matrix of shape = (self.na, self.na * self.n_s). One may figure out the connections via sc_index.

See also

orbital_COHP

orbital resolved COHP analysis of the Green function

orbital_ACOHP

orbital resolved COHP analysis of the spectral function

atom_COHP

atomic COHP analysis of the Green function

atom_COOP(E, kavg=True, isc=None, orbitals=None, uc=False)[source]

Atomic COOP curve of the Green function

Parameters
  • E (float or int) – the energy or the energy index of COOP. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned COOP is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned COOP from unit-cell ([None, None, None]) to the given supercell, the default is all COOP for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain COOP matrix elements for a subset of orbitals, all other are set to 0.

  • uc (bool, optional) – whether the returned COOP are only in the unit-cell. If True this will return a sparse matrix of shape = (self.na, self.na), else, it will return a sparse matrix of shape = (self.na, self.na * self.n_s). One may figure out the connections via sc_index.

See also

orbital_COOP

orbital resolved COOP analysis of the Green function

atom_COOP_from_orbital

transfer an orbital COOP to atomic COOP

atom_ACOOP

atomic COOP analysis of the spectral function

atom_COHP

atomic COHP analysis of the Green function

atom_COOP_from_orbital(COOP, uc=False)[source]

Calculate the atomic COOP curve from the orbital COOP

The atomic COOP are a sum over all orbital COOP:

\[\mathrm{COOP}_{\alpha\beta} = \sum_{\nu\in\alpha}\sum_{\mu\in\beta} \mathrm{COOP}_{\nu\mu}\]
Parameters
  • COOP (scipy.sparse.csr_matrix) – the orbital COOP as retrieved from orbital_COOP or orbital_ACOOP

  • uc (bool, optional) – whether the returned COOP are only in the unit-cell. If True this will return a sparse matrix of shape = (self.na, self.na), else, it will return a sparse matrix of shape = (self.na, self.na * self.n_s). One may figure out the connections via sc_index.

See also

orbital_COOP

orbital resolved COOP analysis of the Green function

orbital_ACOOP

orbital resolved COOP analysis of the spectral function

atom_COOP

atomic COOP analysis of the Green function

atom_current(elec, E, kavg=True, activity=True, orbitals=None)[source]

Atomic current of atoms

Short hand function for calling orbital_current and atom_current_from_orbital.

Parameters
  • elec (str, int) – the electrode of originating electrons

  • E (float or int) – the energy or energy index of the atom current.

  • kavg (bool, int, optional) – whether the returned atomic current is k-averaged, or an explicit (unweighed) k-point is returned

  • activity (bool, optional) – whether the activity current is returned, see atom_current_from_orbital for details.

  • orbitals (array-like or dict, optional) – only retain orbital currents for a subset of orbitals before calculating bond-current Passed directly to orbital_current.

See also

orbital_current

the orbital current between individual orbitals

bond_current_from_orbital

transfer the orbital current to bond current

bond_current

the bond current (orbital current summed over orbitals)

vector_current

an atomic field current for each atom (Cartesian representation of bond-currents)

atom_current_from_orbital(Jij, activity=True)[source]

Atomic current of atoms by passing the orbital current

The atomic current is a single number specifying a figure of the magnitude current flowing through each atom. It is thus not a quantity that can be related to the physical current flowing in/out of atoms but is merely a number that provides an idea of how much current this atom is redistributing.

The atomic current may have two meanings based on these two equations

\[\begin{split}J_\alpha^{|a|} &=\frac{1}{2} \sum_\beta \Big| \sum_{\nu\in \alpha}\sum_{\mu\in \beta} J_{\nu\mu} \Big| \\ J_\alpha^{|o|} &=\frac{1}{2} \sum_\beta \sum_{\nu\in \alpha}\sum_{\mu\in \beta} \big| J_{\nu\mu} \big|\end{split}\]

If the activity current is requested (activity=True) \(J_\alpha^{\mathcal A} = \sqrt{ J_\alpha^{|a|} J_\alpha^{|o|} }\) is returned.

If activity=False \(J_\alpha^{|a|}\) is returned.

For geometries with all atoms only having 1-orbital, they are equivalent.

Generally the activity current is a more rigorous figure of merit for the current flowing through an atom. More so than than the summed absolute atomic current due to the following reasoning. The activity current is a geometric mean of the absolute bond current and the absolute orbital current. This means that if there is an atom with a large orbital current it will have a larger activity current.

Parameters

Examples

>>> Jij = tbt.orbital_current(0, -1.03) # orbital current @ E = -1 eV originating from electrode ``0``
>>> Ja = tbt.atom_current_from_orbital(Jij)
property base_file

File of the current Sile

bloch(elec)

Bloch-expansion coefficients for an electrode

Parameters

elec (str or int) – bloch expansions of electrode

bond_current(elec, E, kavg=True, isc=None, only='+', orbitals=None, uc=False)[source]

Bond-current between atoms (sum of orbital currents)

Short hand function for calling orbital_current and bond_current_from_orbital.

Parameters
  • elec (str, int) – the electrode of originating electrons

  • E (float or int) – A float for energy in eV, int for explicit energy index Unlike orbital_current this may not be None as the down-scaling of the orbital currents may not be equivalent for all energy points.

  • kavg (bool, int, optional) – whether the returned bond current is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned bond currents from the unit-cell ([None, None, None]) (default) to the given supercell. If [None, None, None] is passed all bond currents are returned.

  • only ({'+', '-', 'all'}) – If “+” is supplied only the positive orbital currents are used, for “-”, only the negative orbital currents are used, else return the sum of both. Please see discussion in orbital_current.

  • orbitals (array-like or dict, optional) – only retain orbital currents for a subset of orbitals before calculating bond-current Passed directly to orbital_current.

  • uc (bool, optional) – whether the returned bond-currents are only in the unit-cell. If True this will return a sparse matrix of shape = (self.na, self.na), else, it will return a sparse matrix of shape = (self.na, self.na * self.n_s). One may figure out the connections via sc_index.

Examples

>>> Jij = tbt.orbital_current(0, -1.0, only='+') # orbital current @ E = -1 eV originating from electrode ``0``
>>> Jab1 = tbt.bond_current_from_orbital(Jij)
>>> Jab2 = tbt.bond_current(0, -1.0)
>>> Jab1 == Jab2
True

See also

orbital_current

the orbital current between individual orbitals

bond_current_from_orbital

transfer the orbital current to bond current

atom_current

the atomic current for each atom (scalar representation of bond-currents)

vector_current

an atomic field current for each atom (Cartesian representation of bond-currents)

bond_current_from_orbital(Jij, only='+', uc=False)[source]

Bond-current between atoms (sum of orbital currents) from an external orbital current

Conversion routine from orbital currents into bond currents.

The bond currents are a sum over all orbital currents:

\[J_{\alpha\beta} = \sum_{\nu\in\alpha}\sum_{\mu\in\beta} J_{\nu\mu}\]

where if

  • only='+': only \(J_{\nu\mu} > 0\) are summed onto the corresponding atom,

  • only='-': only \(J_{\nu\mu} < 0\) are summed onto the corresponding atom,

  • only='all': all \(J_{\nu\mu}\) are summed onto the corresponding atom.

Parameters
  • Jij (scipy.sparse.csr_matrix) – the orbital currents as retrieved from orbital_current

  • only ({'+', '-', 'all'}) – If “+” is supplied only the positive orbital currents are used, for “-”, only the negative orbital currents are used, else return both.

  • uc (bool, optional) – whether the returned bond-currents are only in the unit-cell. If True this will return a sparse matrix of shape = (self.na, self.na), else, it will return a sparse matrix of shape = (self.na, self.na * self.n_s). One may figure out the connections via sc_index.

Examples

>>> Jij = tbt.orbital_current(0, -1.0) # orbital current @ E = -1 eV originating from electrode ``0``
>>> Jab = tbt.bond_current_from_orbital(Jij)
>>> Jab[2,3] # bond current between atom 3 and 4

See also

orbital_current

the orbital current between individual orbitals

bond_current

the bond current (orbital current summed over orbitals)

atom_current

the atomic current for each atom (scalar representation of bond-currents)

vector_current

an atomic field current for each atom (Cartesian representation of bond-currents)

btd(elec=None)

Block-sizes for the BTD method in the device/electrode region

Parameters

elec (str or int, optional) – the BTD block sizes for the device (if none), otherwise the downfolding BTD block sizes for the electrode

property cell

Unit cell in file

chemical_potential(elec)

Return the chemical potential associated with the electrode elec

close()
current(elec_from=0, elec_to=1, kavg=True)[source]

Current from from to to using the k-weights and energy spacings in the file.

Calculates the current as:

\[I(\mu_t - \mu_f) = \frac{e}{h}\int\!\mathrm{d}E\, T(E) [n_F(\mu_t, k_B T_t) - n_F(\mu_f, k_B T_f)]\]

The chemical potential and the temperature are taken from this object.

Parameters
  • elec_from (str, int, optional) – the originating electrode

  • elec_to (str, int, optional) – the absorbing electrode (different from elec_from)

  • kavg (bool, int, optional) – whether the returned current is k-averaged, or an explicit (unweighed) k-point is returned

See also

current_parameter

to explicitly set the electronic temperature and chemical potentials

chemical_potential

routine that defines the chemical potential of the queried electrodes

kT

routine that defines the electronic temperature of the queried electrodes

current_parameter(elec_from, mu_from, kt_from, elec_to, mu_to, kt_to, kavg=True)[source]

Current from from to to using the k-weights and energy spacings in the file.

Calculates the current as:

\[I(\mu_t - \mu_f) = \frac{e}{h}\int\!\mathrm{d}E\, T(E) [n_F(\mu_t, k_B T_t) - n_F(\mu_f, k_B T_f)]\]

The chemical potential and the temperature are passed as arguments to this routine.

Parameters
  • elec_from (str, int) – the originating electrode

  • mu_from (float) – the chemical potential of the electrode (in eV)

  • kt_from (float) – the electronic temperature of the electrode (in eV)

  • elec_to (str, int) – the absorbing electrode (different from elec_from)

  • mu_to (float) – the chemical potential of the electrode (in eV)

  • kt_to (float) – the electronic temperature of the electrode (in eV)

  • kavg (bool, int, optional) – whether the returned current is k-averaged, or an explicit (unweighed) k-point is returned

See also

current

which calculates the current with the chemical potentials and temperatures set in the TBtrans calculation

density_matrix(E, kavg=True, isc=None, orbitals=None, geometry=None)[source]

Density matrix from the Green function at energy E (1/eV)

The density matrix can be used to calculate the LDOS in real-space.

The \(\mathrm{LDOS}(E, \mathbf r)\) may be calculated using the density routine. Basically the LDOS in real-space may be calculated as

\[\rho_{\mathbf G}(E, \mathbf r) = -\frac{1}{\pi}\sum_{\nu\mu}\phi_\nu(\mathbf r)\phi_\mu(\mathbf r) \Im[\mathbf G_{\nu\mu}(E)]\]

where \(\phi\) are the orbitals. Note that the broadening used in the TBtrans calculations ensures the broadening of the density, i.e. it should not be necessary to perform energy averages over the density matrices.

Parameters
  • E (float or int) – the energy or the energy index of density matrix. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned density matrix is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned density matrix from unit-cell ([None, None, None]) to the given supercell, the default is all density matrix elements for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain density matrix elements for a subset of orbitals, all other are set to 0.

  • geometry (Geometry, optional) – geometry that will be associated with the density matrix. By default the geometry contained in this file will be used. However, then the atomic species are probably incorrect, nor will the orbitals contain the basis-set information required to generate the required density in real-space.

See also

Adensity_matrix

spectral function density matrix

Returns

object containing the Geometry and the density matrix elements

Return type

DensityMatrix

dir_file(filename=None, filename_base='')

File of the current Sile

property elecs

List of electrodes

electron_temperature(elec)

Electron bath temperature [Kelvin]

Parameters

elec (str or int) – electrode to extract the temperature from

See also

kT

bath temperature in [eV]

eta(elec=None)

The imaginary part used when calculating the self-energies in eV (or for the device

Parameters

elec (str, int, optional) – electrode to extract the eta value from. If not specified (or None) the device region eta will be returned.

fano(elec_from=0, elec_to=1, kavg=True, zero_T=1e-06)[source]

The Fano-factor for the calculation (requires calculated transmission eigenvalues)

Calculate the Fano factor defined as (or through the shot-noise):

\[\begin{split}F(E) &= \frac{\sum_{k,n} T_{k,n}(E)[1 - T_{k,n}(E)] w_k}{\sum_{k,n} T_{k,n}(E) w_k} \\ &= S(E, V) / S_P(E, V)\end{split}\]

Notes

The default zero_T may change in the future. This calculation will only work for non-polarized calculations since the divisor needs to be the spin-sum. The current implementation uses the full transmission as the divisor.

Examples

For a spin-polarized calculation one should calculate the Fano factor as:

>>> up = get_sile('siesta.TBT_UP.nc')
>>> down = get_sile('siesta.TBT_DN.nc')
>>> fano = up.fano() * up.transmission() + down.fano() * down.transmission()
>>> fano /= up.transmission() + down.transmission()
Parameters
  • elec_from (str, int, optional) – the originating electrode

  • elec_to (str, int, optional) – the absorbing electrode (different from elec_from)

  • kavg (bool, int, optional) – whether the returned Fano factor is k-averaged, or an explicit (unweighed) k-point is returned. In any case the divisor will always be the k-averaged transmission.

  • zero_T (float, optional) – any transmission eigen value lower than this value will be treated as exactly 0.

See also

shot_noise

shot-noise term (zero temperature limit)

noise_power

temperature dependent noise power

property file

File of the current Sile

property geom

The associated geometry from this file

property geometry

The associated geometry from this file

info(elec=None)[source]

Information about the calculated quantities available for extracting in this file

Parameters

elec (str or int) – the electrode to request information from

iter(group=True, dimension=True, variable=True, levels=- 1, root=None)

Iterator on all groups, variables and dimensions.

This iterator iterates through all groups, variables and dimensions in the Dataset

The generator sequence will _always_ be:

  1. Group

  2. Dimensions in group

  3. Variables in group

As the dimensions are generated before the variables it is possible to copy groups, dimensions, and then variables such that one always ensures correct dependencies in the generation of a new SileCDF.

Parameters
  • group (bool (True)) – whether the iterator yields Group instances

  • dimension (bool (True)) – whether the iterator yields Dimension instances

  • variable (bool (True)) – whether the iterator yields Variable instances

  • levels (int (-1)) – number of levels to traverse, with respect to root variable, i.e. number of sub-groups this iterator will return.

  • root (str (None)) – the base root to start iterating from.

Examples

Script for looping and checking each instance.

>>> for gv in self.iter():
...     if self.isGroup(gv):
...         # is group
...     elif self.isDimension(gv):
...         # is dimension
...     elif self.isVariable(gv):
...         # is variable
property k

Sampled k-points in file

kT(elec)

Electron bath temperature [eV]

Parameters

elec (str or int) – electrode to extract the temperature from

See also

electron_temperature

bath temperature in [K]

kindex(k)

Return the index of the k-point that is closests to the queried k-point (in reduced coordinates)

Parameters

k (array_like of float or int) – the queried k-point in reduced coordinates \(]-0.5;0.5]\). If int return it-self.

property kpt

Sampled k-points in file

property lasto

Last orbital of corresponding atom

mu(elec)

Return the chemical potential associated with the electrode elec

property nE

Number of energy-points in file

n_btd(elec=None)

Number of blocks in the BTD partioning

Parameters

elec (str or int, optional) – if None the number of blocks in the device region BTD matrix. Else the number of BTD blocks in the electrode down-folding.

property na

Returns number of atoms in the cell

property na_b

Number of atoms in the buffer region

property na_buffer

Number of atoms in the buffer region

property na_d

Number of atoms in the device region

property na_dev

Number of atoms in the device region

na_down(elec)

Number of atoms in the downfolding region (without device downfolded region)

Parameters

elec (str or int) – Number of downfolding atoms for electrode elec

property na_u

Returns number of atoms in the cell

property ne

Number of energy-points in file

property nk

Number of k-points in file

property nkpt

Number of k-points in file

property no

Returns number of orbitals in the cell

property no_d

Number of orbitals in the device region

no_down(elec)

Number of orbitals in the downfolding region (plus device downfolded region)

Parameters

elec (str or int) – Number of downfolding orbitals for electrode elec

no_e(elec)

Number of orbitals in the downfolded region of the electrode in the device

Parameters

elec (str or int) – Specify the electrode to query number of downfolded orbitals

property no_u

Returns number of orbitals in the cell

noise_power(elec_from=0, elec_to=1, kavg=True)[source]

Noise power from to to using the k-weights and energy spacings in the file (temperature dependent)

Calculates the noise power as

\[S(V) = \frac{2e^2}{h}\sum_k\sum_n \int\mathrm d E \big\{T_{k,n}(E)[f_L(1-f_L)+f_R(1-f_R)] + T_{k,n}(E)[1 - T_{k,n}(E)](f_L - f_R)^2\big\} w_k\]

Where \(f_i\) are the Fermi-Dirac distributions for the electrodes.

Raises

SislInfo – If all of the calculated \(T_{k,n}(E)\) values in the file are above 0.001.

Parameters
  • elec_from (str, int, optional) – the originating electrode

  • elec_to (str, int, optional) – the absorbing electrode (different from elec_from)

  • kavg (bool, int, optional) – whether the returned noise-power is k-averaged, or an explicit (unweighed) k-point is returned

See also

fano

the ratio between the quantum mechanial and the classical shot noise.

shot_noise

shot-noise term (zero temperature limit)

norm(atoms=None, orbitals=None, norm='none')[source]

Normalization factor depending on the input

The normalization can be performed in one of the below methods. In the following \(N\) refers to the normalization constant that is to be used (i.e. the divisor):

'none'

\(N=1\)

'all'

\(N\) equals the number of orbitals in the total device region.

'atom'

\(N\) equals the total number of orbitals in the selected atoms. If orbitals is an argument a conversion of orbitals to the equivalent unique atoms is performed, and subsequently the total number of orbitals on the atoms is used. This makes it possible to compare the fraction of orbital DOS easier.

'orbital'

\(N\) is the sum of selected orbitals, if atoms is specified, this is equivalent to the ‘atom’ option.

Parameters
  • atoms (array_like of int or bool, optional) – only return for a given set of atoms (default to all). NOT allowed with orbitals keyword

  • orbitals (array_like of int or bool, optional) – only return for a given set of orbitals (default to all) NOT allowed with atoms keyword

  • norm ({'none', 'atom', 'orbital', 'all'}) – how the normalization of the summed DOS is performed (see norm routine)

o2p(orbitals, elec=None)

Return the pivoting indices (0-based) for the orbitals, possibly on an electrode

Will warn if an orbital requested is not in the device list of orbitals.

Parameters
  • orbitals (array_like or int) – orbital indices (0-based)

  • elec (str or int or None, optional) – electrode to return pivoting indices of (if None it is the device pivoting indices).

property o_dev

Orbital indices (0-based) of device orbitals (sorted)

See also

pivot

retrieve the device orbitals, non-sorted

orbital_ACOHP(elec, E, kavg=True, isc=None, orbitals=None)[source]

Orbital resolved COHP analysis of the spectral function

This will return a sparse matrix, see scipy.sparse.csr_matrix for details. Each matrix element of the sparse matrix corresponds to the COHP of the underlying geometry.

The COHP analysis can be written as:

\[\mathrm{COHP}^{\mathbf A}_{\nu\mu} = \frac{1}{2\pi} \Re\big[\mathbf A_{\nu\mu} \mathbf H_{\nu\mu} \big]\]
Parameters
  • elec (str or int) – the electrode of the spectral function

  • E (float or int) – the energy or the energy index of COHP. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned COHP is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned COHP from unit-cell ([None, None, None]) to the given supercell, the default is all COHP for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain COHP matrix elements for a subset of orbitals, all other are set to 0.

See also

orbital_COHP

orbital resolved COHP analysis of the Green function

atom_COHP_from_orbital

atomic COHP analysis from an orbital COHP

atom_COHP

atomic COHP analysis of the Green function

atom_ACOHP

atomic COHP analysis of the spectral function

orbital_COOP

orbital resolved COOP analysis of the Green function

atom_COOP_from_orbital

transfer an orbital COOP to atomic COOP

atom_COOP

atomic COOP analysis of the Green function

orbital_ACOOP

orbital resolved COOP analysis of the spectral function

atom_ACOOP

atomic COOP analysis of the spectral function

orbital_ACOOP(elec, E, kavg=True, isc=None, orbitals=None)[source]

Orbital COOP analysis of the spectral function

This will return a sparse matrix, see csr_matrix for details. Each matrix element of the sparse matrix corresponds to the COOP of the underlying geometry.

The COOP analysis can be written as:

\[\mathrm{COOP}^{\mathbf A}_{\nu\mu} = \frac{1}{2\pi} \Re\big[\mathbf A_{\nu\mu} \mathbf S_{\mu\nu} \big]\]

The sum of the COOP DOS is equal to the DOS:

\[\mathrm{ADOS}_{\nu} = \sum_\mu \mathrm{COOP}^{\mathbf A}_{\nu\mu}\]

One can calculate the (diagonal) balanced COOP analysis, see JPCM 15 (2003), 7751-7761 for details. The DBCOOP is given by:

\[\begin{split}D &= \sum_\nu \mathrm{COOP}^{\mathbf A}_{\nu\nu} \\ \mathrm{DBCOOP}^{\mathbf A}_{\nu\mu} &= \mathrm{COOP}^{\mathbf A}_{\nu\mu} / D\end{split}\]

The BCOOP can be looked up in the reference above.

Parameters
  • elec (str or int) – the electrode of the spectral function

  • E (float or int) – the energy or the energy index of COOP. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned COOP is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned COOP from unit-cell ([None, None, None]) to the given supercell, the default is all COOP for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain COOP matrix elements for a subset of orbitals, all other are set to 0.

Examples

>>> ACOOP = tbt.orbital_ACOOP(0, -1.0) # COOP @ E = -1 eV from ``0`` spectral function
>>> ACOOP[10, 11] # COOP value between the 11th and 12th orbital
>>> ACOOP.sum(1).A[tbt.o_dev, 0] == tbt.ADOS(0, sum=False)[tbt.Eindex(-1.0)]
>>> D = ACOOP.diagonal().sum()
>>> ADBCOOP = ACOOP / D

See also

orbital_COOP

orbital resolved COOP analysis of the Green function

atom_COOP_from_orbital

transfer an orbital COOP to atomic COOP

atom_COOP

atomic COOP analysis of the Green function

atom_ACOOP

atomic COOP analysis of the spectral function

orbital_COHP

orbital resolved COHP analysis of the Green function

atom_COHP_from_orbital

atomic COHP analysis from an orbital COHP

atom_COHP

atomic COHP analysis of the Green function

orbital_ACOHP

orbital resolved COHP analysis of the spectral function

atom_ACOHP

atomic COHP analysis of the spectral function

orbital_COHP(E, kavg=True, isc=None, orbitals=None)[source]

Orbital resolved COHP analysis of the Green function

This will return a sparse matrix, see scipy.sparse.csr_matrix for details. Each matrix element of the sparse matrix corresponds to the COHP of the underlying geometry.

The COHP analysis can be written as:

\[\mathrm{COHP}^{\mathbf G}_{\nu\mu} = \frac{-1}{2\pi} \Im\big[(\mathbf G - \mathbf G^\dagger)_{\nu\mu} \mathbf H_{\mu\nu} \big]\]
Parameters
  • E (float or int) – the energy or the energy index of COHP. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned COHP is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned COHP from unit-cell ([None, None, None]) to the given supercell, the default is all COHP for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain COHP matrix elements for a subset of orbitals, all other are set to 0.

Examples

>>> COHP = tbt.orbital_COHP(-1.0) # COHP @ E = -1 eV
>>> COHP[10, 11] # COHP value between the 11th and 12th orbital

See also

atom_COHP_from_orbital

atomic COHP analysis from an orbital COHP

atom_COHP

atomic COHP analysis of the Green function

orbital_ACOHP

orbital resolved COHP analysis of the spectral function

atom_ACOHP

atomic COHP analysis of the spectral function

orbital_COOP

orbital resolved COOP analysis of the Green function

atom_COOP_from_orbital

transfer an orbital COOP to atomic COOP

atom_COOP

atomic COOP analysis of the Green function

orbital_ACOOP

orbital resolved COOP analysis of the spectral function

atom_ACOOP

atomic COOP analysis of the spectral function

orbital_COOP(E, kavg=True, isc=None, orbitals=None)[source]

Orbital COOP analysis of the Green function

This will return a sparse matrix, see scipy.sparse.csr_matrix for details. Each matrix element of the sparse matrix corresponds to the COOP of the underlying geometry.

The COOP analysis can be written as:

\[\mathrm{COOP}^{\mathbf G}_{\nu\mu} = \frac{-1}{2\pi} \Im\big[(\mathbf G - \mathbf G^\dagger)_{\nu\mu} \mathbf S_{\mu\nu} \big]\]

The sum of the COOP DOS is equal to the DOS:

\[\mathrm{DOS}_{\nu} = \sum_\mu \mathrm{COOP}^{\mathbf G}_{\nu\mu}\]

One can calculate the (diagonal) balanced COOP analysis, see JPCM 15 (2003), 7751-7761 for details. The DBCOOP is given by:

\[\begin{split}D &= \sum_\nu \mathrm{COOP}^{\mathbf G}_{\nu\nu} \\ \mathrm{DBCOOP}^{\mathbf G}_{\nu\mu} &= \mathrm{COOP}^{\mathbf G}_{\nu\mu} / D\end{split}\]

The BCOOP can be looked up in the reference above.

Parameters
  • E (float or int) – the energy or the energy index of COOP. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned COOP is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned COOP from unit-cell ([None, None, None]) to the given supercell, the default is all COOP for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • orbitals (array-like or dict, optional) – only retain COOP matrix elements for a subset of orbitals, all other are set to 0.

Examples

>>> COOP = tbt.orbital_COOP(-1.0) # COOP @ E = -1 eV
>>> COOP[10, 11] # COOP value between the 11th and 12th orbital
>>> COOP.sum(1).A[tbt.o_dev, 0] == tbt.DOS(sum=False)[tbt.Eindex(-1.0)]
>>> D = COOP.diagonal().sum()
>>> DBCOOP = COOP / D

See also

atom_COOP_from_orbital

transfer an orbital COOP to atomic COOP

atom_COOP

atomic COOP analysis of the Green function

orbital_ACOOP

orbital resolved COOP analysis of the spectral function

atom_ACOOP

atomic COOP analysis of the spectral function

orbital_COHP

orbital resolved COHP analysis of the Green function

atom_COHP_from_orbital

atomic COHP analysis from an orbital COHP

atom_COHP

atomic COHP analysis of the Green function

orbital_ACOHP

orbital resolved COHP analysis of the spectral function

atom_ACOHP

atomic COHP analysis of the spectral function

orbital_current(elec, E, kavg=True, isc=None, only='all', orbitals=None)[source]

Orbital current originating from elec as a sparse matrix

This will return a sparse matrix, see scipy.sparse.csr_matrix for details. Each matrix element of the sparse matrix corresponds to the orbital indices of the underlying geometry.

When requesting orbital-currents it is vital to consider how the data needs to be analysed before extracting the data. For instance, if only local currents are interesting one should use only='+'. While if one is interested in the transmission between subset of orbitals, only='all' is the correct method.

For inexperienced users it is adviced to try out all three values of only to ensure the correct physics is obtained.

This becomes even more important when the orbital currents are calculated with magnetic fields. With \(\mathbf B\) fields local current loops may form and current does not necessarily flow along the transport direction.

Parameters
  • elec (str, int) – the electrode of originating electrons

  • E (float or int) – the energy or the energy index of the orbital current. If an integer is passed it is the index, otherwise the index corresponding to Eindex(E) is used.

  • kavg (bool, int, optional) – whether the returned orbital current is k-averaged, or an explicit (unweighed) k-point is returned

  • isc (array_like, optional) – the returned bond currents from the unit-cell ([None, None, None]) to the given supercell, the default is all orbital currents for the supercell. To only get unit cell orbital currents, pass [0, 0, 0].

  • only ({'all', '+', '-'}) – which orbital currents to return, all, positive or negative values only. Default to 'all' because it can then be used in the subsequent default arguments for bond_current_from_orbital and atom_current_from_orbital.

  • orbitals (array-like or dict, optional) – only retain orbital currents for a subset of orbitals.

Examples

>>> Jij = tbt.orbital_current(0, -1.0) # orbital current @ E = -1 eV originating from electrode ``0``
>>> Jij[10, 11] # orbital current from the 11th to the 12th orbital
>>> Jij = tbt.orbital_current(0, -1.0,
...     orbitals={tbt.geometry.atoms[0]: [0, 1]})

only retain currents from 1st and 2nd orbitals on first atom type (all atoms of that type in the entire structure.

See also

bond_current_from_orbital

transfer the orbital current to bond current

bond_current

the bond current (orbital current summed over orbitals)

atom_current_from_orbital

transfer the orbital current to atomic current

atom_current

the atomic current for each atom (scalar representation of bond-currents)

vector_current

an atomic field current for each atom (Cartesian representation of bond-currents)

pivot(elec=None, in_device=False, sort=False)

Return the pivoting indices for a specific electrode (in the device region) or the device

Parameters
  • elec (str or int) – the corresponding electrode to return the pivoting indices from

  • in_device (bool, optional) – If True the pivoting table will be translated to the device region orbitals. If sort is also true, this would correspond to the orbitals directly translated to the geometry self.geometry.sub(self.a_dev).

  • sort (bool, optional) – Whether the returned indices are sorted. Mostly useful if you want to handle the device in a non-pivoted order.

Examples

>>> se = tbtncSileTBtrans(...)
>>> se.pivot()
[3, 4, 6, 5, 2]
>>> se.pivot(sort=True)
[2, 3, 4, 5, 6]
>>> se.pivot(0)
[2, 3]
>>> se.pivot(0, in_device=True)
[4, 0]
>>> se.pivot(0, in_device=True, sort=True)
[0, 1]
>>> se.pivot(0, sort=True)
[2, 3]

See also

pivot_down

for the pivot table for electrodes down-folding regions

pivot_down(elec)

Pivoting orbitals for the downfolding region of a given electrode

Parameters

elec (str or int) – the corresponding electrode to get the pivoting indices for

read(*args, **kwargs)

Generic read method which should be overloaded in child-classes

Parameters

kwargs – keyword arguments will try and search for the attribute read_<> and call it with the remaining **kwargs as arguments.

read_data(*args, **kwargs)[source]

Read specific type of data.

This is a generic routine for reading different parts of the data-file.

Parameters
  • geom (bool, optional) – return the geometry

  • atom_current (bool, optional) – return the atomic current flowing through an atom (the activity current)

  • vector_current (bool, optional) – return the orbital currents as vectors

read_geometry(*args, **kwargs)

Returns Geometry object from this file

read_supercell()

Returns SuperCell object from this file

reflection(elec=0, kavg=True, from_single=False)[source]

Reflection into electrode elec

The reflection into electrode elec is calculated as:

\[R(E) = T_{\mathrm{bulk}}(E) - \sum_{\mathrm{to}} T_{\mathrm{elec}\to\mathrm{to}}(E)\]

Another way of calculating the reflection is via:

\[R(E) = T_{\mathrm{bulk}}(E) - \big\{i \mathrm{Tr}[(\mathbf G-\mathbf G^\dagger)\boldsymbol\Gamma_{\mathrm{elec}}] - \mathrm{Tr}[\mathbf G\boldsymbol\Gamma_{\mathrm{elec}}\mathbf G^\dagger\boldsymbol\Gamma_{\mathrm{elec}}]\big\}\]

Both are identical, however, numerically they may be different. Particularly when the bulk transmission is very large compared to the transmission to the other electrodes one should prefer the first equation.

Parameters
  • elec (str, int, optional) – the backscattered electrode

  • kavg (bool, int, optional) – whether the returned reflection is k-averaged, or an explicit (unweighed) k-point is returned

  • from_single (bool, optional) – whether the reflection is calculated using the Green function and a single scattering matrix Eq. (2) above (true), otherwise Eq. (1) will be used (false).

See also

transmission

the total transmission

transmission_eig

the transmission decomposed in eigenchannels

transmission_bulk

the total transmission in a periodic lead

shot_noise(elec_from=0, elec_to=1, classical=False, kavg=True)[source]

Shot-noise term from to to using the k-weights

Calculates the shot-noise term according to classical (also known as the Poisson value). If classical is True the shot-noise calculated is:

\[S_P(E, V) = \frac{2e^2}{h}|V|\sum_k\sum_n T_{k,n}(E) w_k = \frac{2e^3}{h}|V|T(E)\]

while for classical False (default) the Fermi-Dirac statistics is taken into account:

\[S(E, V) = \frac{2e^2}{h}|V|\sum_k\sum_n T_{k,n}(E) [1 - T_{k,n}(E)] w_k\]
Raises

SislInfo – If all of the calculated \(T_{k,n}(E)\) values in the file are above 0.001.

Parameters
  • elec_from (str, int, optional) – the originating electrode

  • elec_to (str, int, optional) – the absorbing electrode (different from elec_from)

  • classical (bool, optional) – which shot-noise to calculate, default to non-classical

  • kavg (bool, int, optional) – whether the returned shot-noise is k-averaged, or an explicit (unweighed) k-point is returned

See also

fano

the ratio between the quantum mechanial and the classical shot noise.

noise_power

temperature dependent noise power

transmission(elec_from=0, elec_to=1, kavg=True)[source]

Transmission from elec_from to elec_to.

The transmission between two electrodes may be retrieved from the Sile.

The transmission is calculated as:

\[T(E) = \mathrm{Tr}[\mathbf{G}\boldsymbol\Gamma_{\mathrm{from}}\mathbf{G}^\dagger\boldsymbol\Gamma_{\mathrm{to}}]\]

where all quantities are energy dependent.

Parameters
  • elec_from (str, int, optional) – the originating electrode

  • elec_to (str, int, optional) – the absorbing electrode (different from elec_from)

  • kavg (bool, int, optional) – whether the returned transmission is k-averaged, or an explicit (unweighed) k-point is returned

See also

transmission_eig

the transmission decomposed in eigenchannels

transmission_bulk

the total transmission in a periodic lead

reflection

total reflection back into the electrode

transmission_bulk(elec=0, kavg=True)[source]

Bulk transmission for the elec electrode

The bulk transmission is equivalent to creating a 2 terminal device with electrode elec tiled 3 times.

Parameters
  • elec (str, int, optional) – the bulk electrode

  • kavg (bool, int, optional) – whether the returned transmission are k-averaged, or an explicit (unweighed) k-point is returned

See also

transmission

the total transmission

transmission_eig

the transmission decomposed in eigenchannels

reflection

total reflection back into the electrode

transmission_eig(elec_from=0, elec_to=1, kavg=True)[source]

Transmission eigenvalues from elec_from to elec_to.

Parameters
  • elec_from (str, int, optional) – the originating electrode

  • elec_to (str, int, optional) – the absorbing electrode (different from elec_from)

  • kavg (bool, int, optional) – whether the returned transmission eigenvalues are k-averaged, or an explicit (unweighed) k-point is returned

See also

transmission

the total transmission

transmission_bulk

the total transmission in a periodic lead

vector_current(elec, E, kavg=True, only='+', orbitals=None)[source]

Vector for each atom describing the mean path for the current travelling through the atom

See vector_current_from_bond for details.

Parameters
  • elec (str or int) – the electrode of originating electrons

  • E (float or int) – the energy or energy index of the vector current. Unlike orbital_current this may not be None as the down-scaling of the orbital currents may not be equivalent for all energy points.

  • kavg (bool, int, optional) – whether the returned vector current is k-averaged, or an explicit (unweighed) k-point is returned

  • only ({'+', '-', 'all'}) – By default only sum outgoing vector currents ('+'). The incoming vector currents may be retrieved by '-', while the average incoming and outgoing direction can be obtained with 'all'. In the last case the vector currents are divided by 2 to ensure the length of the vector is compatibile with the other options given a pristine system.

  • orbitals (array-like or dict, optional) – only retain orbital currents for a subset of orbitals before calculating bond-current Passed directly to orbital_current.

Returns

array of vectors per atom in the Geometry (only non-zero for device atoms)

Return type

numpy.ndarray

See also

orbital_current

the orbital current between individual orbitals

bond_current_from_orbital

transfer the orbital current to bond current

bond_current

the bond current (orbital current summed over orbitals)

atom_current

the atomic current for each atom (scalar representation of bond-currents)

vector_current_from_bond(Jab)[source]

Vector for each atom being the sum of bond-current times the normalized bond between the atoms

The vector current is defined as:

\[\mathbf J_\alpha = \sum_\beta \frac{r_\beta - r_\alpha}{|r_\beta - r_\alpha|} \cdot J_{\alpha\beta}\]

Where \(J_{\alpha\beta}\) is the bond current between atom \(\alpha\) and \(\beta\) and \(r_\alpha\) are the atomic coordinates.

Parameters

Jab (scipy.sparse.csr_matrix) – the bond currents as retrieved from bond_current

Returns

array of vectors per atom in the Geometry (only non-zero for device atoms)

Return type

numpy.ndarray

See also

orbital_current

the orbital current between individual orbitals

bond_current_from_orbital

transfer the orbital current to bond current

bond_current

the bond current (orbital current summed over orbitals)

atom_current

the atomic current for each atom (scalar representation of bond-currents)

property wk

Weights of k-points in file

property wkpt

Weights of k-points in file

write(*args, **kwargs)

Generic write method which should be overloaded in child-classes

Parameters

**kwargs – keyword arguments will try and search for the attribute write_ and call it with the remaining **kwargs as arguments.

write_tbtav(*args, **kwargs)[source]

Convert this to a TBT.AV.nc file, i.e. all k dependent quantites are averaged out.

This command will overwrite any previous file with the ending TBT.AV.nc and thus will not take notice of any older files.

Parameters

file (str) – output filename

property xa

Atomic coordinates in file

property xyz

Atomic coordinates in file