sisl.io.tbtrans.tbtprojncSileTBtrans
- class sisl.io.tbtrans.tbtprojncSileTBtrans(filename, mode='r', lvl=0, access=1, *args, **kwargs)
Bases:
tbtncSileTBtrans
TBtrans projection file object
Plotting
Plotting functions for the
tbtprojncSileTBtrans
class.plot.geometry
(*args[, ...])Calls
read_geometry
and creates aGeometryPlot
from its output.plot.pdos
([geometry, ...])Creates a
PDOSData
object and then plots aPdosPlot
from it.Methods
ADOS
(elec_mol_proj[, E, kavg, atoms, ...])Projected spectral density of states (DOS) (1/eV)
Adensity_matrix
(elec_mol_proj, E[, kavg, ...])Projected 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[, method])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
(E[, elec, kavg, isc, orbitals, uc])Atomic COHP curve of the spectral function
atom_ACOOP
(E[, elec, kavg, isc, orbitals, uc])Atomic COOP curve of the spectral function
atom_current
([elec, elec_other, activity, ...])Atomic current of atoms, a scalar quantity quantifying how much currents flows through an atom
atom_transmission
(E[, elec, activity, kavg, ...])Atomic transmission at energy
E
of atoms, a scalar quantity quantifying how much transmission flows through an atombase_directory
([relative_to])Retrieve the base directory of the file, relative to the path relative_to
bloch
(elec)Bloch-expansion coefficients for an electrode
bond_current
([elec, elec_other, kavg, isc, ...])Bond current between atoms (sum of orbital currents)
bond_transmission
(E[, elec, kavg, isc, ...])Bond transmission between atoms at a specific energy
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
()dir_file
([filename, filename_base])File of the current Sile
eigenstate
(molecule[, k, all])Return the eigenstate on the projected molecule
electron_temperature
(elec)Electron bath temperature [Kelvin]
eta
([elec])The imaginary part used when calculating the self-energies in eV (or for the device
info
([molecule])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
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_mol_proj, E[, kavg, isc, ...])Orbital COHP analysis of the projected spectral function
orbital_ACOOP
(elec_mol_proj, E[, kavg, isc, ...])Orbital COOP analysis of the projected spectral function
orbital_current
([elec, elec_other, kavg, ...])Orbital current originating from elec as a sparse matrix
orbital_transmission
(E, elec_mol_proj, ...)Transmission at energy
E
between orbitals originating from elecpivot
([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
projections
(molecule)List of projections on molecule
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
Returns Lattice object from this file
reflection
([elec, kavg, from_single])Reflection into electrode elec
Reduce an atomic sparse matrix to a vector contribution of each atom
sparse_orbital_to_atom
(Dij[, uc, sum_dup])Reduce a sparse matrix in orbital sparse to a sparse matrix in atomic indices
sparse_orbital_to_scalar
(Dij[, activity])Atomic scalar contribution of atoms for a sparse orbital matrix
sparse_orbital_to_vector
(Dij[, uc, sum_dup])Reduce an orbital sparse matrix to a vector contribution of each atom
transmission
(elec_mol_proj_from, ...[, kavg])Transmission from mol_proj_elec_from to mol_proj_elec_to
transmission_bulk
([elec, kavg])Bulk transmission for the elec electrode
transmission_eig
(elec_mol_proj_from, ...[, kavg])Transmission eigenvalues from elec_mol_proj_from to elec_mol_proj_to
vector_current
([elec, elec_other, kavg, ...])Vector for each atom being the sum of bond currents times the normalized bond vector between the atoms
vector_transmission
(E[, elec, kavg, isc, ...])Vector for each atom being the sum of bond transmissions times the normalized bond vector 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.
Attributes
Sampled energy-points in file
Atomic indices (0-based) of device atoms
Atomic indices (0-based) of device atoms (sorted)
File of the current Sile
Unit cell in file
List of electrodes
File of the current Sile
Same as
geometry
, but deprecatedThe associated geometry from this file
Sampled k-points in file
Sampled k-points in file
Last orbital of corresponding atom
List of regions where state projections may happen
Number of energy-points in file
Returns number of atoms in the cell
Number of atoms in the buffer region
Number of atoms in the buffer region
Number of atoms in the device region
Number of atoms in the device region
Returns number of atoms in the cell
Number of energy-points in file
Number of k-points in file
Number of k-points in file
Returns number of orbitals in the cell
Number of orbitals in the device region
Returns number of orbitals in the cell
Orbital indices (0-based) of device orbitals (sorted)
Weights of k-points in file
Weights of k-points in file
Atomic coordinates in file
Atomic coordinates in file
- ADOS(elec_mol_proj, E: EType | None = None, kavg=True, atoms=None, orbitals=None, sum: bool = True, norm: NormType = 'none')[source]
Projected spectral density of states (DOS) (1/eV)
Extract the projected 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_{i\in \{I\}} [\mathbf{G}(E)|i\rangle\langle i|\Gamma_\mathfrak{el}|i\rangle\langle i|\mathbf{G}^\dagger]_{ii}(E)\]where \(|i\rangle\) may be a sum of states. The normalization constant (\(N\)) is defined in the routine
norm
and depends on the arguments.- Parameters:
elec_mol_proj (
str
ortuple
) – originating projected spectral function (<electrode>.<molecule>.<projection>)E (
float
orint
, optional) – optionally only return the DOS of atoms at a given energy pointkavg (
bool
,int
orarray_like
, optional) – whether the returned DOS is k-averaged, an explicit k-point or a selection of k-pointsatoms (
array_like
ofint
orbool
, optional) – only return for a given set of atoms (default to all). NOT allowed with orbital keywordorbitals (
array_like
ofint
orbool
, optional) – only return for a given set of orbitals (default to all) NOT allowed with atoms keywordsum – whether the returned quantities are summed or returned as is, i.e. resolved per atom/orbital.
norm – how the normalization of the summed DOS is performed (see
norm
routine).
- Adensity_matrix(elec_mol_proj, E: EType, kavg=True, isc=None, orbitals=None, geometry: Geometry | None = None)[source]
Projected spectral function density matrix at energy
E
(1/eV)The projected 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\[\boldsymbol\rho_{\mathbf A_{\mathfrak{el}}}(E, \mathbf r) = \frac{1}{2\pi}\sum_{ij}\phi_i(\mathbf r)\phi_j(\mathbf r) \Re[\mathbf A_{\mathfrak{el}, ij}(E)]\]where \(\phi\) are the real-space 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_mol_proj (
str
ortuple
) – the projected electrode of originating electronsE – the density matrix corresponding to the energy.
kavg (
bool
,int
orarray_like
, optional) – whether the returned density matrix is k-averaged, an explicit k-point or a selection of k-pointsisc (
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
ordict
, optional) – only retain density matrix elements for a subset of orbitals, all other are set to 0.geometry – 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.
- Returns:
DensityMatrix (
the object containing the Geometry
andthe density matrix elements
)
- BDOS(elec: str | int = 0, E: str | float | None = None, kavg=True, sum: bool = True, norm: Literal['none', 'atom', 'orbital', 'all'] = 'none') ndarray
Bulk density of states (DOS) (1/eV).
Extract the bulk DOS from electrode elec.
\[\mathrm{BDOS}_\mathfrak{el}(E) = -\frac{1}{\pi} \Im\mathbf{G}(E)\]This returns the density of states for the full (Bloch-expanded) electrode. When
norm
is ‘none’, the DOS is the full DOS for all electrode atoms (fully expanded), if you want to get the DOS for the minimal (un-expanded) electrode unit-cell, then divide bynp.prod(tbt.bloch(elec))
. Whennorm
is anything else, it will be normalised to the number of atoms/orbitals in the electrode.- Parameters:
elec – electrode where the bulk DOS is returned
E – 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 returnedsum – whether the returned quantities are summed or returned as is, i.e. resolved per atom/orbital.
norm – whether the returned quantities are summed over all orbitals or normed by number of orbitals in the electrode. Currently one cannot extract DOS per atom/orbital.
- DOS(E: str | float | None = None, kavg=True, atoms=None, orbitals=None, sum: bool = True, norm: Literal['none', 'atom', 'orbital', 'all'] = 'none') ndarray
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_{i\in \{I\}} \Im \mathbf{G}_{ii}(E)\]The normalization constant (\(N\)) is defined in the routine
norm
and depends on the arguments.- Parameters:
E – 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 returnedatoms (
array_like
ofint
orbool
, 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
ofint
orbool
, 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 – whether the returned quantities are summed or returned as is, i.e. resolved per atom/orbital.
norm – how the normalization of the summed DOS is performed (see
norm
routine)
- Eindex(E: Etype, method: Literal['nearest', 'above', 'below'] = 'nearest')
Return the closest energy index corresponding to the energy
E
- Parameters:
E – 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.
method – how non-equal values should be located. * nearest takes the closests value * above takes the closests value above
E
. * below takes the closests value belowE
.
- 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
orint
) – atomic indices (0-based)
- a_down(elec: str | int, bulk: bool = False)
Down-folding atomic indices for a given electrode
- Parameters:
elec – electrode to retrieve indices for
bulk – whether the returned indices are only in the pristine electrode, or the down-folding region (electrode + downfolding region, not in device)
- a_elec(elec: str | int)
Electrode atomic indices for the full geometry (sorted)
- Parameters:
elec – electrode to retrieve indices for
- atom_ACOHP(E: str | float, elec: str | int = 0, kavg=True, isc=None, orbitals=None, uc: bool = False) csr_matrix
Atomic COHP curve of the spectral function
- Parameters:
E – the atomic COHP corresponding to the energy.
elec – the electrode of the spectral function
kavg (
bool
,int
, optional) – whether the returned COHP is k-averaged, or an explicit (unweighed) k-point is returnedisc (
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
ordict
, optional) – only retain COHP matrix elements for a subset of orbitals, all other are set to 0.uc – whether the returned COHP are only in the unit-cell. If
True
this will return a sparse matrix ofshape = (self.na, self.na)
, else, it will return a sparse matrix ofshape = (self.na, self.na * self.n_s)
. One may figure out the connections viasc_index
.
See also
orbital_COOP
orbital resolved COOP analysis of the Green function
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
atomic COHP analysis of the Green function
orbital_ACOHP
orbital resolved COHP analysis of the spectral function
- atom_ACOOP(E: str | float, elec: str | int = 0, kavg=True, isc=None, orbitals=None, uc: bool = False) csr_matrix
Atomic COOP curve of the spectral function
The atomic COOP are a sum over all orbital COOP:
\[\mathrm{COOP}_{IJ} = \sum_{i\in I}\sum_{j\in J} \mathrm{COOP}_{ij}\]This is a shorthand for calling
orbital_ACOOP
andsparse_orbital_to_atom
in order.- Parameters:
E – the atomic COOP corresponding to the energy.
elec – the electrode of the spectral function
kavg (
bool
,int
, optional) – whether the returned COOP is k-averaged, or an explicit (unweighed) k-point is returnedisc (
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
ordict
, optional) – only retain COOP matrix elements for a subset of orbitals, all other are set to 0.uc – whether the returned COOP are only in the unit-cell. If
True
this will return a sparse matrix ofshape = (self.na, self.na)
, else, it will return a sparse matrix ofshape = (self.na, self.na * self.n_s)
. One may figure out the connections viasc_index
.
See also
orbital_COOP
orbital resolved COOP analysis of the Green function
atom_COOP
atomic COOP analysis of the Green function
orbital_ACOOP
orbital resolved COOP analysis of the spectral function
orbital_COHP
orbital resolved COHP analysis of the Green function
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
- atom_current(elec: str | int = 0, elec_other: str | int = 1, activity: bool = True, kavg=True, isc=None, orbitals=None) ndarray
Atomic current of atoms, a scalar quantity quantifying how much currents flows through an atom
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}\mathbf j_I^{|a|} &=\frac 12 \sum_{\{J\}} \Big| \sum_{i\in I}\sum_{j\in J} \mathbf J_{ij} \Big| \\ \mathbf j_I^{|o|} &=\frac 12 \sum_{i\in I}\sum_{j\in\{J\}} \big| J_{ij} \big|\end{split}\]\[\]If the activity is requested (
activity=True
) \(\mathbf j_I^{\mathcal A} = \sqrt{\mathbf j_I^{|a|} \mathbf j_I^{|o|} }\) is returned.If
activity=False
\(\mathbf j_I^{|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:
elec – the originating electrode
elec_other – this electrode determines the other chemical potential. As such the orbital currents does not reflect the current going from elec to elec_other!
activity –
True
to return the activity current, see explanation abovekavg (
bool
,int
, optional) – whether the returned orbital current is k-averaged, or an explicit (unweighed) k-point is returnedisc (
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]
.orbitals (
array-like
ordict
, optional) – only retain orbital currents for a subset of orbitals.
Examples
>>> Jij = tbt.orbital_current(0, 1, what="all") # orbital current originating from electrode ``0`` >>> Ja = tbt.sparse_orbital_to_scalar(Jij)
Notes
Calculating the current between two electrodes with the same chemical potential will return a matrix filled with 0’s since there is no bias window.
The currents does not reflect the current going from elec_from to elec_other!
See also
orbital_transmission
energy resolved transmission between orbitals
orbital_current
bias window integrated transmissions
bond_transmission
energy resolved transmissions between atoms
bond_current
bias window integrated transmissions (orbital current summed over orbitals)
vector_transmission
an atomic field transmission for each atom (Cartesian representation of bond-transmissions)
vector_current
an atomic field current for each atom (Cartesian representation of bond-currents)
atom_transmission
energy resolved atomic transmission for each atom (scalar representation of bond-transmissions)
- atom_transmission(E: str | float, elec: str | int = 0, activity: bool = True, kavg=True, isc=None, orbitals=None) ndarray
Atomic transmission at energy
E
of atoms, a scalar quantity quantifying how much transmission flows through an atomThe atomic transmission is a single number specifying a figure of the magnitude transmission flowing through each atom. It is thus not a quantity that can be related to the physical transmission flowing in/out of atoms but is merely a number that provides an idea of how much this atom is redistributing.
The atomic transmission may have two meanings based on these two equations
\[\begin{split}T_I^{|a|} &=\frac 12 \sum_{\{J\}} \Big| \sum_{i\in I}\sum_{j\in J} \mathbf T_{ij} \Big| \\ T_I^{|o|} &=\frac 12 \sum_{i\in I}\sum_{j\in\{J\}} \big| T_{ij} \big|\end{split}\]\[\]If the activity is requested (
activity=True
) \(T_I^{\mathcal A} = \sqrt{T_I^{|a|} T_I^{|o|} }\) is returned. If the activity current is requested (activity=True
)If
activity=False
\(T_I^{|a|}\) is returned.For geometries with all atoms only having 1-orbital, they are equivalent.
Generally the activity is a more rigorous figure of merit for the transmission flowing through an atom. More so than than the summed absolute atomic transmission due to the following reasoning. The activity transmission is a geometric mean of the absolute bond transmission and the absolute orbital transmission. This means that if there is an atom with a large orbital transmission it will have a larger activity.
- Parameters:
E – the atomic transmission corresponding to the energy.
elec – the originating electrode
activity –
True
to return the activity, see explanation abovekavg (
bool
,int
, optional) – whether the returned atomic transmissions are k-averaged, or an explicit (unweighed) k-point is returnedisc (
array_like
, optional) – the returned transmissions from the unit-cell ([None, None, None]
) to the given supercell, the default is all orbital transmissions are used for the supercell. To only get unit cell orbital currents, pass[0, 0, 0]
.orbitals (
array-like
ordict
, optional) – only retain orbital currents for a subset of orbitals.
Examples
>>> Jij = tbt.orbital_transmission(-1., what"all") # transmission @ E=-1 eV from electrode ``0`` >>> Ja = tbt.sparse_orbital_to_scalar(Jij)
See also
orbital_transmission
energy resolved transmission between orbitals
orbital_current
bias window integrated transmissions
bond_transmission
energy resolved transmissions between atoms
bond_current
bias window integrated transmissions (orbital current summed over orbitals)
vector_transmission
an atomic field transmission for each atom (Cartesian representation of bond-transmissions)
vector_current
an atomic field current for each atom (Cartesian representation of bond-currents)
atom_current
the atomic current for each atom (scalar representation of bond-currents)
- base_directory(relative_to='.')
Retrieve the base directory of the file, relative to the path relative_to
- bloch(elec: str | int)
Bloch-expansion coefficients for an electrode
- Parameters:
elec – bloch expansions of electrode
- bond_current(elec: str | int = 0, elec_other: str | int = 1, kavg=True, isc=None, what: str = 'all', orbitals=None, uc: bool = False) csr_matrix
Bond current between atoms (sum of orbital currents)
Short hand function for calling
orbital_current
andsparse_orbital_to_atom
.The bond currents are a sum over all orbital currents:
\[J_{IJ} = \sum_{i\in I}\sum_{j\in J} J_{ij}\]- Parameters:
elec – the electrode of originating electrons
elec_other – this electrode determines the other chemical potential. As such the orbital currents does not reflect the current going from elec to elec_other!
kavg (
bool
,int
, optional) – whether the returned bond current is k-averaged, or an explicit (unweighed) k-point is returnedisc (
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.what (
{"all"/"both"/"+-"/"inout", "+"/"out", "-"/"in"}
) – If +/out is supplied only the positive currents are used (going out) for -/in, only the negative currents are used (going in), else return both. Please see discussion inorbital_current
.orbitals (
array-like
ordict
, optional) – only retain currents for a subset of orbitals before calculating bond current Passed directly toorbital_current
.uc – whether the returned currents are only in the unit-cell (supercell currents will be folded to their unit-cell equivalents). If True this will return a sparse matrix of
shape = (self.na, self.na)
, else, it will return a sparse matrix ofshape = (self.na, self.na * self.n_s)
. One may figure out the connections viasc_index
.
Examples
>>> Jij = tbt.orbital_current(0, 1, what="out") # orbital current originating from electrode ``0`` >>> Jab1 = tbt.sparse_orbital_to_atom(Jij) >>> Jab2 = tbt.bond_current(0, 1, what="out") >>> Jab1 == Jab2 True
Notes
Calculating the current between two electrodes with the same chemical potential will return a matrix filled with 0’s since there is no bias window.
The currents does not reflect the current going from elec_from to elec_other!
See also
orbital_transmission
energy resolved transmission between orbitals
orbital_current
bias window integrated transmissions
bond_transmission
energy resolved transmissions between atoms
vector_transmission
an atomic field transmission for each atom (Cartesian representation of bond-transmissions)
vector_current
an atomic field current for each atom (Cartesian representation of bond-currents)
atom_transmission
energy resolved atomic transmission for each atom (scalar representation of bond-transmissions)
atom_current
the atomic current for each atom (scalar representation of bond-currents)
- bond_transmission(E: str | float, elec: str | int = 0, kavg=True, isc=None, what: str = 'all', orbitals=None, uc: bool = False) csr_matrix
Bond transmission between atoms at a specific energy
Short hand function for calling
orbital_transmission
andsparse_orbital_to_atom
.The bond transmissions are a sum over all orbital transmissions
\[T_{IJ}(E) = \sum_{i\in I}\sum_{j\in J} T_{ij}(E)\]- Parameters:
E – the bond transmission corresponding to the energy.
elec – the electrode of originating electrons
kavg (
bool
,int
, optional) – whether the returned bond transmissions is k-averaged, or an explicit (unweighed) k-point is returnedisc (
array_like
, optional) – the returned transmissions from the unit-cell ([None, None, None]
) (default) to the given supercell. If[None, None, None]
is passed all transmissions are returned.what (
{"all"/"both"/"+-"/"inout", "+"/"out", "-"/"in"}
) – If +/out is supplied only the positive transmissions are used (going out) for -/in, only the negative transmissions are used (going in), else return both. Please see discussion inorbital_transmission
.orbitals (
array-like
ordict
, optional) – only retain transmissions for a subset of orbitals before calculating bond transmissions Passed directly toorbital_transmission
.uc – whether the returned transmissions are only in the unit-cell (supercell bonds will be folded to their unit-cell equivalents). If True this will return a sparse matrix of
shape = (self.na, self.na)
, else, it will return a sparse matrix ofshape = (self.na, self.na * self.n_s)
. One may figure out the connections viasc_index
.
Examples
>>> Jij = tbt.orbital_transmission(-1.0, what="out") # orbital transmission @ E = -1 eV originating from electrode ``0`` >>> Jab1 = tbt.sparse_orbital_to_atom(Jij)[ >>> Jab2 = tbt.bond_transmission(-1.0, what="out") >>> Jab1 == Jab2 True
See also
orbital_transmission
energy resolved transmission between orbitals
orbital_current
bias window integrated transmissions
bond_current
bias window integrated transmissions (orbital current summed over orbitals)
vector_transmission
an atomic field transmission for each atom (Cartesian representation of bond-transmissions)
vector_current
an atomic field current for each atom (Cartesian representation of bond-currents)
atom_transmission
energy resolved atomic transmission for each atom (scalar representation of bond-transmissions)
atom_current
the atomic current for each atom (scalar representation of bond-currents)
- btd(elec: int | str | None = None)
Block-sizes for the BTD method in the device/electrode region
- Parameters:
elec – the BTD block sizes for the device (if none), otherwise the downfolding BTD block sizes for the electrode
- chemical_potential(elec: str | int) float
Return the chemical potential associated with the electrode elec
- close()
- dir_file(filename=None, filename_base='')
File of the current Sile
- eigenstate(molecule, k=None, all=True)[source]
Return the eigenstate on the projected molecule
The eigenstate object will contain the geometry as the parent object. The eigenstate will be in the Lowdin basis:
\[|\psi'_i\rangle = \mathbf S^{1/2} |\psi_i\rangle\]
- electron_temperature(elec: str | int) float
Electron bath temperature [Kelvin]
- Parameters:
elec – electrode to extract the temperature from
See also
kT
bath temperature in [eV]
- eta(elec: int | str | None = None) float
The imaginary part used when calculating the self-energies in eV (or for the device
- Parameters:
elec – electrode to extract the eta value from. If not specified (or None) the device region eta will be returned.
- info(molecule=None)[source]
Information about the calculated quantities available for extracting in this file
- 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:
Group
Dimensions in group
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 instancesdimension (
bool
(True)) – whether the iterator yields Dimension instancesvariable (
bool
(True)) – whether the iterator yields Variable instanceslevels (
int
(-1)) – number of levels to traverse, with respect toroot
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
- kT(elec: str | int) float
Electron bath temperature [eV]
- Parameters:
elec – 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)
- n_btd(elec: int | str | None = None) int
Number of blocks in the BTD partioning
- Parameters:
elec – if None the number of blocks in the device region BTD matrix. Else the number of BTD blocks in the electrode down-folding.
- na_down(elec: str | int) int
Number of atoms in the downfolding region (without device downfolded region)
- Parameters:
elec – Number of downfolding atoms for electrode elec
- no_down(elec: str | int) int
Number of orbitals in the downfolding region (plus device downfolded region)
- Parameters:
elec – Number of downfolding orbitals for electrode elec
- no_e(elec: str | int) int
Number of orbitals in the downfolded region of the electrode in the device
- Parameters:
elec – Specify the electrode to query number of downfolded orbitals
- norm(atoms=None, orbitals=None, norm: Literal['none', 'atom', 'orbital', 'all'] = 'none') int
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
ofint
orbool
, optional) – only return for a given set of atoms (default to all). NOT allowed with orbitals keywordorbitals (
array_like
ofint
orbool
, optional) – only return for a given set of orbitals (default to all) NOT allowed with atoms keywordnorm – how the normalization of the summed DOS is performed (see
norm
routine)
- o2p(orbitals, elec: int | str | None = 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
orint
) – orbital indices (0-based)elec – electrode to return pivoting indices of (if None it is the device pivoting indices).
- orbital_ACOHP(elec_mol_proj, E: str | float, kavg=True, isc=None, orbitals=None)[source]
Orbital COHP analysis of the projected 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}_{ij} = \frac{1}{2\pi} \Re\big[\mathbf A_{ij} \mathbf H_{ij} \big]\]- Parameters:
elec_mol_proj (
str
ortuple
) – the electrode of the projected spectral functionE – the COHP values corresponding to the energy.
kavg (
bool
,int
orarray_like
, optional) – whether the returned COHP is k-averaged, an explicit k-point or a selection of k-pointsisc (
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
ordict
, optional) – only retain COHP matrix elements for a subset of orbitals, all other are set to 0.
See also
atom_COHP_from_orbital
atomic COHP analysis from an orbital COHP
atom_ACOHP
atomic COHP analysis of the projected spectral function
atom_COOP_from_orbital
transfer an orbital COOP to atomic COOP
orbital_ACOOP
orbital resolved COOP analysis of the projected spectral function
atom_ACOOP
atomic COOP analysis of the projected spectral function
- orbital_ACOOP(elec_mol_proj, E: str | float, kavg=True, isc=None, orbitals=None)[source]
Orbital COOP analysis of the projected 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 COOP of the underlying geometry.The COOP analysis can be written as:
\[\mathrm{COOP}^{\mathbf A}_{ij} = \frac{1}{2\pi} \Re\big[\mathbf A_{ij} \mathbf S_{ji} \big]\]The sum of the COOP DOS is equal to the DOS:
\[\mathrm{ADOS}_{i} = \sum_j \mathrm{COOP}^{\mathbf A}_{ij}\]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_i \mathrm{COOP}^{\mathbf A}_{ii} \\ \mathrm{DBCOOP}^{\mathbf A}_{ij} &= \mathrm{COOP}^{\mathbf A}_{ij} / D\end{split}\]The BCOOP can be looked up in the reference above.
- Parameters:
elec_mol_proj (
str
ortuple
) – the electrode of the spectral functionE – the COOP values corresponding to the energy.
kavg (
bool
,int
orarray_like
, optional) – whether the returned COOP is k-averaged, an explicit k-point or a selection of k-pointsisc (
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
ordict
, optional) – only retain COOP matrix elements for a subset of orbitals, all other are set to 0.
Examples
>>> ACOOP = tbt.orbital_ACOOP('Left.C60.HOMO', -1.0) # COOP @ E = -1 eV from ``Left.C60.HOMO`` 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
atom_COOP_from_orbital
transfer an orbital COOP to atomic COOP
atom_ACOOP
atomic COOP analysis of the projected spectral function
atom_COHP_from_orbital
atomic COHP analysis from an orbital COHP
orbital_ACOHP
orbital resolved COHP analysis of the projected spectral function
atom_ACOHP
atomic COHP analysis of the projected spectral function
- orbital_current(elec: str | int = 0, elec_other: str | int = 1, kavg=True, isc=None, what: str = 'all', orbitals=None) csr_matrix
Orbital current originating from elec as a sparse matrix
This is the bias window integrated quantity of
orbital_transmission
. As such it represents how the current is flowing at an applied bias from a given electrode.\[J_{ij} = \frac eh\int_{\mu_1}^{\mu_2} \!\mathrm dE\, T_{ij} [n_F(\mu_2, k_B T_2) - n_F(\mu_1, k_B T_1)]\]with \(T_{\langle\rangle}\) being the electronic temperature of the respective reservoir.
- Parameters:
elec – the originating electrode
elec_other – this electrode determines the other chemical potential. As such the orbital currents does not reflect the current going from elec to elec_other!
kavg (
bool
,int
, optional) – whether the returned orbital current is k-averaged, or an explicit (unweighed) k-point is returnedisc (
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]
.what (
{"all"/"both"/"+-"/"inout", "+"/"out", "-"/"in"}
) – which orbital currents to return, all, positive (outgoing) or negative (incoming). Default to"all"
because it can then be used in the subsequent default arguments forsparse_orbital_to_atom
andsparse_orbital_to_scalar
.orbitals (
array-like
ordict
, optional) – only retain orbital currents for a subset of orbitals.
Notes
Calculating the current between two electrodes with the same chemical potential will return a matrix filled with 0’s since there is no bias window.
The currents does not reflect the current going from elec_from to elec_other!
See also
orbital_transmission
energy resolved transmission between orbitals
bond_transmission
energy resolved transmissions between atoms
bond_current
bias window integrated transmissions (orbital current summed over orbitals)
vector_transmission
an atomic field transmission for each atom (Cartesian representation of bond-transmissions)
vector_current
an atomic field current for each atom (Cartesian representation of bond-currents)
atom_transmission
energy resolved atomic transmission for each atom (scalar representation of bond-transmissions)
atom_current
the atomic current for each atom (scalar representation of bond-currents)
- orbital_transmission(E: str | float, elec_mol_proj, *args, **kwargs)[source]
Transmission at energy
E
between orbitals originating from elecEach matrix element of the sparse matrix corresponds to the orbital indices of the underlying geometry (including buffer and electrode atoms).
When requesting orbital-transmissions it is vital to consider how the data needs to be analysed before extracting the data. For instance, if only local transmission pathways are interesting one should use
what="+"
to retain the positive orbital transmissions. While if one is interested in the transmission between subset of orbitals,what="all"
is the correct method to account for loop transmissions.The orbital transmissions are calculated as described in the TBtrans manual:
\[T_{ij}(E) = i [ (\mathbf H_{ji} - E\mathbf S_{ji}) \mathbf A_{ij}(E) - (\mathbf H_{ij} - E\mathbf S_{ij}) \mathbf A_{ji}(E)],\]It is easy to show that the above matrix obeys \(T_{ij}=-T_{ji}\).
For inexperienced users it is adviced to try out all three values of
what
to ensure the correct physics is obtained.This becomes even more important when the orbital transmissions are calculated with magnetic fields. With \(\mathbf B\) fields local transmission loops may form and the pathways does not necessarily flow along the transport direction.
For correct interpretation of the orbital transmissions it is vital that one integrates the full Brillouin zone without any symmetry operations, see Section 5.4 in [10].
- Parameters:
E – the orbital transmission corresponding to the energy.
elec – the electrode of originating electrons
kavg (
bool
,int
, optional) – whether the returned orbital transmission is k-averaged, or an explicit (unweighed) k-point is returnedisc (
array_like
, optional) – the returned transmissions from the unit-cell ([None, None, None]
) to the given supercell, the default is all transmissions for the supercell. To only get unit cell transmissions, pass[0, 0, 0]
.what (
{"all"/"both"/"+-"/"inout", "+"/"out", "-"/"in"}
) – which transmissions to return, all, positive (outgoing) or negative (incoming).orbitals (
array-like
ordict
, optional) – only retain transmissions for a subset of orbitals (including their supercell equivalents)
- Returns:
A `scipy.sparse.csr_matrix
containing the supercell transmission pathways`,or
orbital transmissions.
Examples
>>> Jij = tbt.orbital_transmission(-1.0) # orbital current @ E = -1 eV originating from electrode ``0`` >>> Jij[10, 11] # orbital transmission from the 11th to the 12th orbital
>>> Jij = tbt.orbital_transmission(-1.0, ... orbitals={tbt.geometry.atoms[0]: [0, 1]})
only retain transmissions from 1st and 2nd orbitals on first atom type (all atoms of that type in the entire structure.
See also
orbital_current
bias window integrated transmissions
bond_transmission
energy resolved transmissions between atoms
bond_current
bias window integrated transmissions (orbital current summed over orbitals)
vector_transmission
an atomic field transmission for each atom (Cartesian representation of bond-transmissions)
vector_current
an atomic field current for each atom (Cartesian representation of bond-currents)
atom_transmission
energy resolved atomic transmission for each atom (scalar representation of bond-transmissions)
atom_current
the atomic current for each atom (scalar representation of bond-currents)
- pivot(elec: int | str | None = None, in_device: bool = False, sort: bool = False)
Return the pivoting indices for a specific electrode (in the device region) or the device
- Parameters:
elec – Can be None, to specify the device region pivot indices (default). Otherwise, it corresponds to the pivoting indicies in the downfolding region.
in_device – 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 geometryself.geometry.sub(self.a_dev)
.sort – 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: str | int)
Pivoting orbitals for the downfolding region of a given electrode
- Parameters:
elec – the corresponding electrode to get the pivoting indices for
- plot.geometry(*args, data_kwargs={}, axes: Axes = ['x', 'y', 'z'], atoms: AtomsIndex = None, atoms_style: Sequence[AtomsStyleSpec] = [], atoms_scale: float = 1.0, atoms_colorscale: Colorscale | None = None, drawing_mode: Literal['scatter', 'balls', None] = None, bind_bonds_to_ats: bool = True, points_per_bond: int = 20, bonds_style: StyleSpec = {}, bonds_scale: float = 1.0, bonds_colorscale: Colorscale | None = None, show_atoms: bool = True, show_bonds: bool = True, show_cell: Literal['box', 'axes', False] = 'box', cell_style: StyleSpec = {}, nsc: tuple[int, int, int] = (1, 1, 1), atoms_ndim_scale: tuple[float, float, float] = (16, 16, 1), bonds_ndim_scale: tuple[float, float, float] = (1, 1, 10), dataaxis_1d: np.ndarray | Callable | None = None, arrows: Sequence[AtomArrowSpec] = (), backend='plotly') GeometryPlot
Calls
read_geometry
and creates aGeometryPlot
from its output.- Parameters:
axes – The axes to project the geometry to.
atoms – The atoms to plot. If None, all atoms are plotted.
atoms_style – List of style specifications for the atoms. See the showcase notebooks for examples.
atoms_scale – Scaling factor for the size of all atoms.
atoms_colorscale – Colorscale to use for the atoms in case the color attribute is an array of values. If None, the default colorscale is used for each backend.
drawing_mode – The method used to draw the atoms.
bind_bonds_to_ats – Whether to display only bonds between atoms that are being displayed.
points_per_bond – When the points are drawn using points instead of lines (e.g. in some frameworks to draw multicolor bonds), the number of points used per bond.
bonds_style – Style specification for the bonds. See the showcase notebooks for examples.
bonds_scale – Scaling factor for the width of all bonds.
bonds_colorscale – Colorscale to use for the bonds in case the color attribute is an array of values. If None, the default colorscale is used for each backend.
show_atoms – Whether to display the atoms.
show_bonds – Whether to display the bonds.
show_cell – Mode to display the cell. If False, the cell is not displayed.
cell_style – Style specification for the cell. See the showcase notebooks for examples.
nsc – Number of unit cells to display in each direction.
atoms_ndim_scale – Scaling factor for the size of the atoms for different dimensionalities (1D, 2D, 3D).
bonds_ndim_scale – Scaling factor for the width of the bonds for different dimensionalities (1D, 2D, 3D).
dataaxis_1d – Only meaningful for 1D plots. The data to plot on the Y axis.
arrows – List of arrow specifications to display. See the showcase notebooks for examples.
backend – The backend to use to generate the figure.
See also
GeometryPlot
The plot class used to generate the plot.
read_geometry
The method called to get the data.
- plot.pdos(geometry: Geometry | None = None, *, groups: Sequence[OrbitalStyleQuery] = [{'name': 'DOS'}], Erange: tuple[float, float] = (-2, 2), E_axis: Literal['x', 'y'] = 'x', line_mode: Literal['line', 'scatter', 'area_line'] = 'line', line_scale: float = 1.0, backend: str = 'plotly') PdosPlot
Creates a
PDOSData
object and then plots aPdosPlot
from it.- Parameters:
geometry – Full geometry of the system (including scattering and electrode regions). Right now only used to get the basis of each atom, which is not stored in the TBT.nc file.
groups – List of orbital specifications to filter and accumulate the PDOS. The contribution of each group will be displayed in a different line. See showcase notebook for examples.
Erange – The energy range to plot.
E_axis – Axis to project the energies.
line_mode – Mode used to draw the PDOS lines.
line_scale – Scaling factor for the width of all lines.
backend – The backend to generate the figure.
See also
PdosPlot
The plot class used to generate the plot.
PDOSData
The class to which data is converted.
- projections(molecule)[source]
List of projections on molecule
- Parameters:
molecule (
str
) – name of molecule to retrieve projections on
- 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)
Read specific type of data.
This is a generic routine for reading different parts of the data-file.
- Parameters:
geometry (
bool
, optional) – return the geometryvector_transmission (
bool
, optional) – return the bond transmissions as vectorsvector_current (
bool
, optional) – return the bond currents as vectorsatom_transmission (
bool
, optional) – return the atomic transmission flowing through an atom (the activity current)atom_current (
bool
, optional) – return the atomic current flowing through an atom (the activity current)
- read_geometry(*args, **kwargs)
Returns Geometry object from this file
- reflection(elec: str | int = 0, kavg=True, from_single: bool = False) ndarray
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 – the backscattered electrode
kavg (
bool
,int
, optional) – whether the returned reflection is k-averaged, or an explicit (unweighed) k-point is returnedfrom_single – 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
- sparse_atom_to_vector(Dab) ndarray
Reduce an atomic sparse matrix to a vector contribution of each atom
Notes
This routine may be moved to a
sisl.utility
at some point since it would be a generic routine usable for other parts ofsisl
.- Parameters:
Dab (
scipy.sparse.csr_matrix
) – the input sparse matrix in atomic indices
- sparse_orbital_to_atom(Dij, uc: bool = False, sum_dup: bool = True) csr_matrix
Reduce a sparse matrix in orbital sparse to a sparse matrix in atomic indices
This algorithm may keep the same non-zero entries, but will return a new csr_matrix with duplicate indices.
Notes
This routine may be moved to a
sisl.utility
at some point since it would be a generic routine usable for other parts ofsisl
.- Parameters:
Dij (
scipy.sparse.csr_matrix
) – the input sparse matrix in orbital formatuc – whether the returned data are only in the unit-cell. If
True
this will return a sparse matrix ofshape = (self.na, self.na)
, else, it will return a sparse matrix ofshape = (self.na, self.na * self.n_s)
. One may figure out the connections viasc_index
.sum_dup – duplicates will be summed if this is true, in this case, no duplicates are present in the returned sparse matrix. If false, duplicates may exist for multi-orbital systems.
- sparse_orbital_to_scalar(Dij, activity: bool = True) ndarray
Atomic scalar contribution of atoms for a sparse orbital matrix
The atomic contribution is a single number specifying a figure of the magnitude of sparse matrix elements for each atom. It is thus not a quantity that can be related to any physical quantity that the sparse matrix may represent but is merely a number that provides an idea of how much this atom is governing the data in the matrix.
The atomic contribution may have two meanings based on these two equations
\[\begin{split}\mathbf a_I^{|a|} &=\frac 12 \sum_{\{J\}} \Big| \sum_{i\in I}\sum_{j\in J} \mathbf A_{ij} \Big| \\ \mathbf a_I^{|o|} &=\frac 12 \sum_{i\in I}\sum_{j\in\{J\}} \big| A_{ij} \big|\end{split}\]If the activity is requested (
activity=True
) \(\mathbf a_I^{\mathcal A} = \sqrt{\mathbf a_I^{|a|} \mathbf a_I^{|o|} }\) is returned.If
activity=False
\(\mathbf a_I^{|a|}\) is returned.For geometries with all atoms only having 1-orbital, they are equivalent.
- Parameters:
Dij (
scipy.sparse.csr_matrix
) – the orbital sparse matrix.activity –
True
to return the atomic activity, see explanation above
Notes
This routine may be moved to a
sisl.utility
at some point since it would be a generic routine usable for other parts ofsisl
.Examples
>>> Jij = tbt.orbital_current(0, -1.03, what="both") # orbital current @ E = -1 eV originating from electrode ``0`` >>> Ja = tbt.sparse_orbital_to_scalar(Jij)
- sparse_orbital_to_vector(Dij, uc: bool = False, sum_dup: bool = True) ndarray
Reduce an orbital sparse matrix to a vector contribution of each atom
Equivalent to calling
sparse_orbital_to_atom
andsparse_atom_to_vector
.Notes
This routine may be moved to a
sisl.utility
at some point since it would be a generic routine usable for other parts ofsisl
.- Parameters:
Dij (
scipy.sparse.csr_matrix
) – the input sparse matrixuc – whether the returned data are only in the unit-cell. If
True
this will return a sparse matrix ofshape = (self.na, self.na)
, else, it will return a sparse matrix ofshape = (self.na, self.na * self.n_s)
. One may figure out the connections viasc_index
.sum_dup – duplicates will be summed if this is true, in this case, no duplicates are present in the returned sparse matrix. If false, duplicates may exist for multi-orbital systems.
- transmission(elec_mol_proj_from, elec_mol_proj_to, kavg=True)[source]
Transmission from mol_proj_elec_from to mol_proj_elec_to
- Parameters:
elec_mol_proj_from (
str
ortuple
) – the originating scattering projection (<electrode>.<molecule>.<projection>)elec_mol_proj_to (
str
ortuple
) – the absorbing scattering projection (<electrode>.<molecule>.<projection>)kavg (
bool
,int
orarray_like
, optional) – whether the returned transmission is k-averaged, an explicit k-point or a selection of k-points
See also
transmission_eig
projected transmission decomposed in eigenchannels
- transmission_bulk(elec: str | int = 0, kavg=True) ndarray
Bulk transmission for the elec electrode
The bulk transmission is equivalent to creating a 2 terminal device with electrode elec tiled 3 times.
- Parameters:
See also
transmission
the total transmission
transmission_eig
the transmission decomposed in eigenchannels
reflection
total reflection back into the electrode
- transmission_eig(elec_mol_proj_from, elec_mol_proj_to, kavg=True)[source]
Transmission eigenvalues from elec_mol_proj_from to elec_mol_proj_to
- Parameters:
elec_mol_proj_from (
str
ortuple
) – the originating scattering projection (<electrode>.<molecule>.<projection>)elec_mol_proj_to (
str
ortuple
) – the absorbing scattering projection (<electrode>.<molecule>.<projection>)kavg (
bool
,int
orarray_like
, optional) – whether the returned transmission is k-averaged, an explicit k-point or a selection of k-points
See also
transmission
projected transmission
- vector_current(elec: str | int = 0, elec_other: str | int = 1, kavg=True, isc=None, what: str = 'all', orbitals=None) ndarray
Vector for each atom being the sum of bond currents times the normalized bond vector between the atoms
The vector current is defined as:
\[\mathbf J_I = \sum_J \frac{\mathbf r^{(J)} - \mathbf r^{(I)}}{|\mathbf r^{(J)} - \mathbf r^{(I)}|} \cdot J_{IJ}\]Where \(J_{IJ}\) is the bond current between atom \(I\) and \(J\) and \(\mathbf r^{(\langle\rangle)}\) are the atomic coordinates.
- Parameters:
elec – the electrode of originating electrons
elec_other – this electrode determines the other chemical potential. As such the vector currents does not reflect the current going from elec to elec_other!
kavg (
bool
,int
, optional) – whether the returned vector current is k-averaged, or an explicit (unweighed) k-point is returnedisc (
array_like
, optional) – the returned currents from the unit-cell ([None, None, None]
) to the given supercell, the default is all currents for the supercell. To only get unit cell orbital currents, pass[0, 0, 0]
.what (
{"all"/"both"/"+-"/"inout", "+"/"out", "-"/"in"}
) – The outgoing currents may be retrieved by"out"
. The incoming currents may be retrieved by"in"
, while the average incoming and outgoing direction can be obtained with"both"
. In the last case the vector currents are divided by 2 to ensure the length of the vector is compatible with the other options given a pristine system.orbitals (
array-like
ordict
, optional) – only retain currents for a subset of orbitals before calculating currents Passed directly toorbital_current
.
Notes
Calculating the current between two electrodes with the same chemical potential will return a matrix filled with 0’s since there is no bias window.
The currents does not reflect the current going from elec_from to elec_other!
- Returns:
numpy.ndarray
– array of vectors per atom in the Geometry (only non-zero for device atoms)
See also
orbital_transmission
energy resolved transmission between orbitals
orbital_current
bias window integrated transmissions
bond_transmission
energy resolved transmissions between atoms
bond_current
bias window integrated transmissions (orbital current summed over orbitals)
vector_transmission
an atomic field transmission for each atom (Cartesian representation of bond-transmissions)
atom_transmission
energy resolved atomic transmission for each atom (scalar representation of bond-transmissions)
atom_current
the atomic current for each atom (scalar representation of bond-currents)
- vector_transmission(E: str | float, elec: str | int = 0, kavg=True, isc=None, what='all', orbitals=None) ndarray
Vector for each atom being the sum of bond transmissions times the normalized bond vector between the atoms
The vector transmission is defined as:
\[\mathbf T_I = \sum_J \frac{\mathbf r^{(J)} - \mathbf r^{(I)}}{|\mathbf r^{(J)} - \mathbf r^{(I)}|} \cdot T_{IJ}\]Where \(T_{IJ}\) is the bond transmission between atom \(I\) and \(J\) and \(\mathbf r^{(\langle\rangle)}\) are the atomic coordinates.
- Parameters:
E – the vector transmission corresponding to the energy.
elec – the electrode of originating electrons
kavg (
bool
,int
, optional) – whether the returned vector transmission is k-averaged, or an explicit (unweighed) k-point is returnedisc (
array_like
, optional) – the returned vectors from the unit-cell ([None, None, None]
) to the given supercell, the default is all vectors for the supercell. To only get unit cell vectors, pass[0, 0, 0]
.what (
{"all"/"both"/"+-"/"inout", "+"/"out", "-"/"in"}
) – The outgoing vectors may be retrieved by"out"
. The incoming vectors may be retrieved by"in"
, while the average incoming and outgoing direction can be obtained with"both"
. In the last case the vector transmissions are divided by 2 to ensure the length of the vector is compatible with the other options; given a pristine system.orbitals (
array-like
ordict
, optional) – only retain transmissions for a subset of orbitals before calculating bond transmissions Passed directly toorbital_transmission
.
- Returns:
numpy.ndarray
– array of vectors per atom in the Geometry (only non-zero for device atoms)
See also
orbital_transmission
energy resolved transmission between orbitals
orbital_current
bias window integrated transmissions
bond_transmission
energy resolved transmissions between atoms
bond_current
bias window integrated transmissions (orbital current summed over orbitals)
vector_current
an atomic field current for each atom (Cartesian representation of bond-currents)
atom_transmission
energy resolved atomic transmission for each atom (scalar representation of bond-transmissions)
atom_current
the atomic current for each atom (scalar representation of bond-currents)
- 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)
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 a_buf
Atomic indices (0-based) of device atoms
- property a_dev
Atomic indices (0-based) of device atoms (sorted)
- atom_COHP = None
- atom_COOP = None
- property base_file
File of the current Sile
- current = None
- current_parameter = None
- density_matrix = None
- property elecs
List of electrodes
- fano = None
- property file
File of the current Sile
- property molecules
List of regions where state projections may happen
- noise_power = None
- property o_dev
Orbital indices (0-based) of device orbitals (sorted)
See also
pivot
retrieve the device orbitals, non-sorted
- orbital_COHP = None
- orbital_COOP = None
- plot
Plotting functions for the
tbtprojncSileTBtrans
class.
- shot_noise = None