TransmissionSpectrum

class TransmissionSpectrum(configuration, energies=None, kpoints=None, kpoints_weights=None, self_energy_calculator=None, energy_zero_parameter=None, infinitesimal=None, device_transmission_method=None, enforce_zero_in_band_gap=None)

Constructor for the TransmissionSpectrum object.

Parameters:
  • configuration (DeviceConfiguration | BulkConfiguration) – The configuration with an attached calculator for which to calculate the transmission spectrum. A BulkConfiguration must be constructed as an electrode, i.e. its unit cell vectors A and B are orthogonal to C, and C is in the z-direction.

  • energies (PhysicalQuantity of type energy.) – The energies for which to calculate the transmission spectrum.
    Default: Energy range that covers the bias window, plus \(30 k_B T\).

  • kpoints (sequence (size 3) of int | MonkhorstPackGrid | KpointDensity | RegularKpointGrid | AdaptiveGrid | sequence of sequence (size 3) of float) – The k-points for which to calculate the transmission spectrum. This can either be given as a k-point grid or a list of fractional k-points, e.g., [[0.0, 0.0, 0.0], [0.0, 0.0, 0.1], ...]. Note that the k-points must be in the xy-plane.
    Default: MonkhorstPackGrid(na, nb) where (na, nb) is the sampling used for the self-consistent calculation.

  • kpoints_weights (sequence of float) – The weight of each k-point.
    Default: The weights corresponding to the MonkhorstPackGrid, or a list of [1.0, 1.0, …] if k-points are specified as floats.

  • self_energy_calculator (DirectSelfEnergy | RecursionSelfEnergy | SparseRecursionSelfEnergy | KrylovSelfEnergy) – The SelfEnergyCalculator to be used for the transmission spectrum.
    Default: RecursionSelfEnergy(storage_strategy=NoStorage())

  • energy_zero_parameter (AverageFermiLevel | AbsoluteEnergy) – Specifies the choice for the energy zero.
    Default: AverageFermiLevel

  • infinitesimal (PhysicalQuantity of type energy.) – Small positive energy, used to move the transmission calculation away from the real axis. This is only relevant for recursion-style self-energy calculators.
    Default: 1.0e-6*eV

  • device_transmission_method (GreensFunction | SparseGreensFunction.) – The method employed in a device transmission spectrum calculation. This parameter is only used if configuration is a DeviceConfiguration. Note that GreensFunction in TransmissionSpectrum is only supported for a maximum of 2 (1) processes_per_contour_point for Unpolarized, Polarized (Noncollinear, SpinOrbit) calculations.
    Default: GreensFunction(processes_per_contour_point=1).

  • enforce_zero_in_band_gap (bool) – Flag which specifies whether the transmission values are enforced to zero inside the band gap. This parameter is only used if configuration is a DeviceConfiguration.
    Default: True

adaptive()
Returns:

True if adaptive grids are used, otherwise False.

Return type:

bool

adaptiveGrid(spin_index, energy_index, tolerance=None, error_measure=None, maximum_number_of_levels=None)

Method for construction of a new AdaptiveGrid, which can be used for restarting from a previous step.

Parameters:
  • spin_index (Non-negative int.) – The spin index for which to extract the transmission values, k-points and triangles for the new adaptive grid.
    Default: 0

  • energy_index (Non-negative int.) – The energy index for which to extract the transmission values, k-points and triangles for the new adaptive grid.
    Default: 0

  • tolerance (Positive float.) – Tolerance for adaptive grid integration.
    Default: 1e-6

  • error_measure (Absolute | Relative) – Method to calculate error in adaptive grid integration. Either Absolute error) or Relative error.
    Default: Absolute

  • maximum_number_of_levels (Non-negative int.) – Maximum number of refinement steps in the adaptive grid.
    Default: 10

Returns:

The adaptive grid object constructed with the transmission values, k-points and triangles corresponding to the spin and energy index.

Return type:

AdaptiveGrid

bias(positive_current_convention=None)
Parameters:

positive_current_convention (LeftToRight | RightToLeft) – The convention for the direction of the positive current. If LeftToRight the bias is defined as V(left) - V(right), whereas for RightToLeft the bias is defined as V(right) - V(left).
Default: LeftToRight

Returns:

The applied bias.

Return type:

PhysicalQuantity of type voltage.

conductance(electrode_voltages=None, electrode_temperatures=None, coupling_rates=None, spin=None)

Calculate the differential conductance.

Parameters:
  • electrode_voltages (PhysicalQuantity of type voltage.) – The left and right electrode voltages to be used for defining the bias window for the current calculation. If no voltage is specified it will use the voltages in the self-consistent calculation. If voltages are specified the current will be obtained non-self-consistent.
    Default: The electrode voltages from the calculator.

  • electrode_temperatures (PhysicalQuantity of type temperature.) – The electrode temperatures to be used in the current calculation.
    Default: The temperatures from the calculator.

  • coupling_rates (list of two floats) – Two fractions (alpha_L, alpha_R) with alpha_L + alpha_R = 1, which specifies the fraction of the differential voltage, alpha_L*dV, added to the left electrode voltage, and alpha_R*dV, added to the right electrode voltage.
    Default: [0.5, 0.5]

  • spin (Spin.Up | Spin.Down | Spin.Sum | Spin.X | Spin.Y | Spin.Z | Spin.UpUp | Spin.DownDown | Spin.RealUpDown | Spin.ImagUpDown | Spin.All) – The spin of the current.
    Default: Spin.Sum

Returns:

The differential conductance for the given spin type. For Spin.All a list of (Spin.Sum, Spin.X, Spin.Y, Spin.Z) will be returned.

Return type:

PhysicalQuantity of type conductance.

current(electrode_voltages=None, electrode_temperatures=None, spin=None, positive_current_convention=None)

Calculate the current with positive direction from left to right.

Parameters:
  • electrode_voltages (PhysicalQuantity of type voltage.) – The left and right electrode voltages to be used for defining the bias window for the current calculation. If no voltage is specified it will use the voltages in the self-consistent calculation. If voltages are specified the current will be obtained non-self-consistent.
    Default: The electrode voltages from the calculator.

  • electrode_temperatures (PhysicalQuantity of type temperature.) – The electrode temperatures to be used in the current calculation.
    Default: The temperatures from the calculator.

  • spin (Spin.Up | Spin.Down | Spin.Sum | Spin.X | Spin.Y | Spin.Z | Spin.UpUp | Spin.DownDown | Spin.RealUpDown | Spin.ImagUpDown | Spin.All) – The spin of the current.
    Default: Spin.Sum

  • positive_current_convention (LeftToRight | RightToLeft) – The convention for the direction of the positive current.
    Default: LeftToRight

Returns:

The current for the given spin type. For Spin.All a list of (Spin.Sum, Spin.X, Spin.Y, Spin.Z) will be returned.

Return type:

PhysicalQuantity of type electric current.

electrodeFermiLevels()
Returns:

The Fermi levels of the left and right electrodes in absolute energies.

Return type:

PhysicalQuantity of type energy.

electrodeFermiTemperatures()
Returns:

The Fermi temperature of the left and right electrodes used in this transmission spectrum.

Return type:

PhysicalQuantity of type temperature.

electrodeVoltages()
Returns:

The left and right electrode voltages used in this transmission spectrum.

Return type:

PhysicalQuantity of type voltage.

energies()
Returns:

The energies used in this transmission spectrum

Return type:

PhysicalQuantity of type energy

energyZero()
Returns:

The energy zero used for the energy scale in this transmission spectrum.

Return type:

PhysicalQuantity of type energy.

energyZeroParameter()
Returns:

The specified choice for the energy zero.

Return type:

AverageFermiLevel | AbsoluteEnergy

enforceZeroInBandGap()
Returns:

Whether the transmission values are enforced to zero inside the band gap.

Return type:

bool

evaluate(kpoints=None, spin=None)

Evaluates the k-point averaged transmission spectrum.

Parameters:
  • kpoints (MonkhorstPackGrid | RegularKpointGrid | list(n) of list(3) of float.) – The set of k-points to average over. For calculations with AdaptiveGrid this parameter cannot be set and the average is over all k-points.
    Default: Average over all k-points.

  • spin (Spin.Up | Spin.Down | Spin.Sum | Spin.X | Spin.Y | Spin.Z | Spin.UpUp | Spin.DownDown | Spin.RealUpDown | Spin.ImagUpDown | Spin.All) – The spin component to query result for.
    Default: Spin.Sum

Returns:

The k-point averaged transmission spectrum for the given spin type. For Spin.All a list of (Spin.Sum, Spin.X, Spin.Y, Spin.Z) will be returned.

Return type:

list(4) of array(n_energies) | array(n_energies)

infinitesimal()
Returns:

The infinitesimal used for calculating the transmission.

Return type:

PhysicalQuantity of type energy

kpoints()
Returns:

The k-points used in this transmission spectrum.

Return type:

array(n_kpoints, 3) of float | dict(n_spins, n_energies)

kpointsWeights()
Returns:

The weights of the k-points used in this transmission spectrum.

Return type:

list(n_kpoints) of float | None

metatext()
Returns:

The metatext of the object or None if no metatext is present.

Return type:

str | None

nlinfo()
Returns:

The transmissions spectrum information.

Return type:

dict

nlprint(stream=None)

Print a string containing an ASCII table useful for plotting the AnalysisSpin object.

Parameters:

stream (python stream) – The stream the table should be written to.
Default: NLPrintLogger()

setMetatext(metatext)

Set a given metatext string on the object.

Parameters:

metatext (str | None) – The metatext string that should be set. A value of “None” can be given to remove the current metatext.

spectralConductivity(electrode_voltages=None, electrode_temperatures=None, coupling_rates=None, spin=None, positive_current_convention=None)

Calculate the spectral conductivity.

Parameters:
  • electrode_voltages (PhysicalQuantity of type voltage.) – The left and right electrode voltages to be used for defining the bias window for the calculation.
    Default: The electrode voltages from the calculator.

  • electrode_temperatures (PhysicalQuantity of type temperature.) – The left and right electrode temperatures to be used in the spectral current calculation.
    Default: The temperatures from the calculator.

  • coupling_rates (list of two non-negative float.) – Two fractions (alpha_L, alpha_R) with alpha_L + alpha_R = 1, which specifies the fraction of the differential voltage, alpha_L*dV, added to the left electrode voltage, and alpha_R*dV, added to the right electrode voltage.
    Default: [0.5, 0.5]

  • spin (Spin.Up | Spin.Down | Spin.Sum | Spin.X | Spin.Y | Spin.Z | Spin.UpUp | Spin.DownDown | Spin.RealUpDown | Spin.ImagUpDown | Spin.All) – The spin of the current.
    Default: Spin.Sum

Returns:

The spectral conductivity for the given spin type. For Spin.All a list of (Spin.Sum, Spin.X, Spin.Y, Spin.Z) will be returned.

Return type:

PhysicalQuantity with the unit e**2/hplanck | list of four PhysicalQuantity with the unit e**2/hplanck

spectralCurrent(electrode_voltages=None, electrode_temperatures=None, spin=None, positive_current_convention=None)

Calculate the spectral current.

Parameters:
  • electrode_voltages (PhysicalQuantity of type voltage.) – The left and right electrode voltages to be used for defining the bias window for the current calculation. If no voltage is specified it will use the voltages in the self-consistent calculation. If voltages are specified the current will be obtained non-self-consistent.
    Default: The electrode voltages from the calculator.

  • electrode_temperatures (PhysicalQuantity of type temperature.) – The electrode temperatures to be used in the current calculation.
    Default: The temperatures from the calculator.

  • spin (Spin.Up | Spin.Down | Spin.Sum | Spin.X | Spin.Y | Spin.Z | Spin.UpUp | Spin.DownDown | Spin.RealUpDown | Spin.ImagUpDown | Spin.All) – The spin of the current.
    Default: Spin.Sum

  • positive_current_convention (LeftToRight | RightToLeft) – The convention for the direction of the positive current.
    Default: LeftToRight

Returns:

The spectral current for the given spin type. For Spin.All a list of (Spin.Sum, Spin.X, Spin.Y, Spin.Z) will be returned.

Return type:

PhysicalQuantity with the unit e/hplanck | list of PhysicalQuantity with the unit e/hplanck

spectralThermalCurrent(electrode_voltages=None, electrode_temperatures=None, spin=None, positive_current_convention=None)

Calculate the spectral Thermal current.

Parameters:
  • electrode_voltages (PhysicalQuantity of type voltage.) – The left and right electrode voltages to be used for defining the bias window for the calculation.
    Default: The electrode voltages from the calculator.

  • electrode_temperatures (PhysicalQuantity of type temperature.) – The left and right electrode temperatures to be used in the spectral current calculation.
    Default: The temperatures from the calculator.

  • spin (Spin.Up | Spin.Down | Spin.Sum | Spin.X | Spin.Y | Spin.Z | Spin.UpUp | Spin.DownDown | Spin.RealUpDown | Spin.ImagUpDown | Spin.All) – The spin of the current.
    Default: Spin.Sum

  • positive_current_convention (LeftToRight | RightToLeft) – The convention for the direction of the current. When RightToLeft the thermal current will from the Right electrode.
    Default: LeftToRight

Returns:

The spectral current for the given spin type. For Spin.All a list of (Spin.Sum, Spin.X, Spin.Y, Spin.Z) will be returned.

Return type:

PhysicalQuantity with the unit e/hplanck | list of PhysicalQuantity with the unit e/hplanck

spins()
Returns:

The spins used in this transmission spectrum.

Return type:

tuple(n_spins) of spin component flags.

thermalConductance(temperature=None, fermi_level_shift=None, spin=None)

Calculate the differential conductance.

Parameters:
  • electrode_voltages (PhysicalQuantity of type voltage.) – The left and right electrode voltages to be used for defining the bias window for the conductance calculation. If no voltage is specified it will use the voltages in the self-consistent calculation. If voltages are specified the current will be obtained non-self-consistent.
    Default: The electrode voltages from the calculator.

  • electrode_temperatures (PhysicalQuantity of type temperature.) – The left and right electrode temperatures to be used for the current calculation.
    Default: The temperatures from the calculator.

  • coupling_rates (list of two floats) – Two fractions (alpha_L, alpha_R) with alpha_L + alpha_R = 1, which specifies the fraction of the differential voltage, alpha_L*dV, added to the left electrode voltage, and alpha_R*dV, added to the right electrode voltage.
    Default: [0.5, 0.5]

  • spin (Spin.Up | Spin.Down | Spin.Sum | Spin.X | Spin.Y | Spin.Z | Spin.UpUp | Spin.DownDown | Spin.RealUpDown | Spin.ImagUpDown | Spin.All) – The spin of the conductance.
    Default: Spin.Sum

Returns:

The differential conductance for the given spin type. For Spin.All a list of (Spin.Sum, Spin.X, Spin.Y, Spin.Z) will be returned.

Return type:

PhysicalQuantity of type conductance.

thermalCurrent(electrode_voltages=None, electrode_temperatures=None, spin=None, positive_current_convention=None)

Calculate the Thermal current.

Parameters:
  • electrode_voltages (PhysicalQuantity of type voltage.) – The left and right electrode voltages to be used for defining the bias window for the calculation.
    Default: The electrode voltages from the calculator.

  • electrode_temperatures (PhysicalQuantity of type temperature.) – The left and right electrode temperatures to be used in the spectral current calculation.
    Default: The temperatures from the calculator.

  • spin (Spin.Up | Spin.Down | Spin.Sum | Spin.X | Spin.Y | Spin.Z | Spin.UpUp | Spin.DownDown | Spin.RealUpDown | Spin.ImagUpDown | Spin.All) – The spin of the current.
    Default: Spin.Sum

  • positive_current_convention (LeftToRight | RightToLeft) – The convention for the direction of the positive current.
    Default: RightToLeft

Returns:

The spectral current for the given spin type. For Spin.All a list of (Spin.Sum, Spin.X, Spin.Y, Spin.Z) will be returned.

Return type:

PhysicalQuantity with the unit e/hplanck | list of PhysicalQuantity with the unit e/hplanck

transmission(spin=None)

Determine the transmission coefficients of the transmission spectrum.

Parameters:

spin (Spin.Up | Spin.Down | Spin.Sum | Spin.X | Spin.Y | Spin.Z | Spin.UpUp | Spin.DownDown | Spin.RealUpDown | Spin.ImagUpDown | Spin.All) – The spin component to get the transmission coefficients for.
Default: Spin.All

Returns:

The transmission coefficients for each energy and k-point for the given spin type. For Spin.All a list of (Spin.Sum, Spin.X, Spin.Y, Spin.Z) will be returned.

Return type:

list(4) of array(n_energies, n_kpoints) | array(n_energies, n_kpoints) | list(4) of dict(n_energies, n_kpoints) | dict(n_energies, n_kpoints)

uniqueString()

Return a unique string representing the state of the object.

Usage Examples

Calculate the TransmissionSpectrum and for each energy print out all k-dependent coefficients:

# Define A,B directions of lattice
vector_a = [5.0, 0.0, 0.0]*Angstrom
vector_b = [0.0, 5.0, 0.0]*Angstrom

# Setup a device configuration
electrode = BulkConfiguration(
    bravais_lattice=UnitCell(vector_a, vector_b,
                             [0.0, 0.0, 9.0]*Angstrom),
    elements=[Lithium, Lithium, Lithium],
    cartesian_coordinates=[[ 2.5,  2.5,  1.5],
                           [ 2.5,  2.5,  4.5],
                           [ 2.5,  2.5,  7.5]]*Angstrom
    )
central_region = BulkConfiguration(
    bravais_lattice= UnitCell(vector_a, vector_b,
                             [0.0, 0.0, 22.0]*Angstrom),
    elements=[Lithium, Lithium, Lithium, Hydrogen, Hydrogen,
              Lithium, Lithium,  Lithium],
    cartesian_coordinates=[[  2.5,   2.5,   1.5],
                           [  2.5,   2.5,   4.5],
                           [  2.5,   2.5,   7.5],
                           [  2.5,   2.5,  10.5],
                           [  2.5,   2.5,  11.5],
                           [  2.5,   2.5,  14.5],
                           [  2.5,   2.5,  17.5],
                           [  2.5,   2.5,  20.5]]*Angstrom
     )
device_configuration = DeviceConfiguration(
    central_region,
    [electrode, electrode]
    )

# Setup calculators
numerical_accuracy_parameters = NumericalAccuracyParameters(
    k_point_sampling=(1, 1, 100))

electrode_calculator = HuckelCalculator(
    numerical_accuracy_parameters=numerical_accuracy_parameters)

calculator = DeviceHuckelCalculator(
    numerical_accuracy_parameters=numerical_accuracy_parameters,
    electrode_calculators=[electrode_calculator, electrode_calculator])

device_configuration.setCalculator(calculator)

# Calculate the transmission coefficient
transmission_spectrum = TransmissionSpectrum(
    configuration=device_configuration,
    energies=numpy.linspace(-0.5,0.5,10)*eV,
    kpoints=MonkhorstPackGrid(3,3,1),
    energy_zero_parameter=AverageFermiLevel,
    infinitesimal=1e-06*Units.eV
    )

# Print out all k-dependent transmission coefficients
data = transmission_spectrum.transmission()
energies = transmission_spectrum.energies()
kpoints = transmission_spectrum.kpoints()
spins = ["Up","Down"]
#print data
for s in range(len(spins)):
    print("Transmission coefficient of spin component ", spins[s])
    for i in range(data[s].shape[0]):
        print('Transmission at energy = %12.6f eV' % (energies[i].inUnitsOf(eV)))
        print('      kx         ky      transmission ')
        for j in range(data[s].shape[1]):
            print('%10.4f %10.4f %16.6e' % \
                  (kpoints[j][0],kpoints[j][1],data[s][i][j]))
        print()

li_h2_li_trans.py

Alternatively, using the MonkhorstPackGrid class, we can also specify the k-points and weights directly. The k-points are specified in units of the reciprocal lattice vectors. Note, that except for the Gamma point, the k-points have weight 2 due to use of time-reversal symmetry, i.e \(T(k)=T(-k)\).

transmission_spectrum = TransmissionSpectrum(
    configuration=device_configuration,
    energies=numpy.linspace(-0.5,0.5,10)*eV,
    kpoints=[[0.0,0.0,0.0],[0.0,1./3.,0.0],[1./3.,-1./3.,0.0],
    [1./3.,0.0,0.0],[1./3.,1./3.,0.0]],
    kpoints_weights = [1.,2.,2.,2.,2.],
    energy_zero_parameter=AverageFermiLevel,
    infinitesimal=1e-06*Units.eV
    )

li_h2_li_trans1.py

Finally, the TransmissionSpectrum can be calculated using an AdaptiveGrid of k-points which is automatically refined in the regions of k-space where a fine k-point sampling is needed:

#----------------------------------------
# Adaptive Grid.
#----------------------------------------

adaptive_grid = AdaptiveGrid(
	kA_range=[-0.5, 0.5],
	kB_range=[-0.5, 0.5],
	tolerance=1e-2,
	error_measure=Relative,
	number_of_initial_levels=2,
	maximum_number_of_levels=7)

#----------------------------------------
# Transmission Spectrum.
#----------------------------------------

transmission_spectrum = TransmissionSpectrum(
    configuration=device_configuration,
    energies=numpy.linspace(-0.5,0.5,10)*eV,
    kpoints=adaptive_grid,
    energy_zero_parameter=AverageFermiLevel,
    infinitesimal=1e-06*Units.eV
    )

trans_adaptive.py

The following lines can be appended to the script to evaluate the current and differential conductance for different electron temperatures in the electrodes:

temp = numpy.linspace(0,1000,21)
currents = []
conductances = []

for t in temp:
    current = transmission_spectrum.current(
                   electrode_temperatures=[200.,t]*Kelvin
		   )
    currents.append(current)

    conductance = transmission_spectrum.conductance(
                       electrode_temperatures=[200.,t]*Kelvin
		       )
    conductances.append(conductance)

print("Temperature[Left]   Temperature[Right]"\
      + "     Current     Differential-Conductance ")
output_format = "      %6.0f K          %6.0f K     %s  %s"
for i in range(len(temp)):
    print(output_format%(200., temp[i], currents[i], conductances[i]))

current.py

A TransmissionSpectrum can also be calculated for a BulkConfiguration. This example show the calculation of the TransmissionSpectrum for a Bulk boron wire:

# -------------------------------------------------------------
# Bulk Configuration
# -------------------------------------------------------------

# Set up lattice
vector_a = [10.0, 0.0, 0.0]*Angstrom
vector_b = [0.0, 10.0, 0.0]*Angstrom
vector_c = [0.0, 0.0, 10.0]*Angstrom
lattice = UnitCell(vector_a, vector_b, vector_c)

# Define elements
elements = [Boron, Boron, Boron, Boron, Boron]

# Define coordinates
fractional_coordinates = [[ 0.5,  0.5,  0. ],
                          [ 0.5,  0.5,  0.2],
                          [ 0.5,  0.5,  0.4],
                          [ 0.5,  0.5,  0.6],
                          [ 0.5,  0.5,  0.8]]

# Set up configuration
bulk_configuration = BulkConfiguration(
    bravais_lattice=lattice,
    elements=elements,
    fractional_coordinates=fractional_coordinates
    )

# -------------------------------------------------------------
# Calculator
# -------------------------------------------------------------
calculator = HuckelCalculator()

bulk_configuration.setCalculator(calculator)
nlprint(bulk_configuration)
bulk_configuration.update()
nlsave('boron_wire_transmission_spectrum.nc', bulk_configuration)

# -------------------------------------------------------------
# Transmission Spectrum
# -------------------------------------------------------------
transmission_spectrum = TransmissionSpectrum(
    configuration=bulk_configuration,
    energies=numpy.linspace(-2,2,101)*eV,
    kpoints=MonkhorstPackGrid(1,1),
    energy_zero_parameter=AverageFermiLevel,
    infinitesimal=1e-06*eV,
    self_energy_calculator=RecursionSelfEnergy(),
    )
nlsave('boron_wire_transmission_spectrum.nc', transmission_spectrum)
nlprint(transmission_spectrum)

boron_wire_transmission_spectrum.py

Notes

To export the data of a TransmissionSpectrum, use the method nlsave.

Unit of transmission

A fully transmitting channel contributes \(1\) per spin.

Note

If you sum the spin channels, the total transmission of an ideal system is 2. Often, it is the energy dependent conductance which is reported in articles,

\[G(E) = \frac{e^2}{h} T(E).\]

If \(G(E)\) is plotted in units of the conductance quantum \(G_0= 2 \frac{e^2}{h}\) for an ideal system,

\[G(E) = 1 G_0.\]

In QuantumATK we usually report the transmission coefficient per spin channel, even for unpolarized systems. Thus, the reported transmission coefficient per spin channel for an ideal system is 1.

See also, TransmissionEigenvalues.

Current

The current is calculated from the transmission coefficient using

\[I(V_L, V_R, T_L, T_R) = \frac{e}{h} \sum_\sigma \int T_\sigma(E) \left[f\left(\frac{E-\mu_R}{k_B T_R}\right) - f\left(\frac{E-\mu_L}{k_B T_L}\right) \right] dE,\]

where \(f\) is the Fermi function, \(T_{L/R}\) is the electron temperatures of the left/right electrode, and \(T_\sigma(E)\) is the transmission coefficient for the spin component \(\sigma\).

The chemical potentials of the left electrode, \(\mu_L = E_F^L-e V_L\), and the right electrode, \(\mu_R = E_F^L-e V_R\), are defined relative to the Fermi level of the left electrode and related to the applied bias through

\[\mu_R - \mu_L = e V_{bias},\]

thus

\[V_{bias} = V_L - V_R .\]

Note

The Fermi levels of the left and right electrodes are not necessarily the same a priori, if for instance the electrodes consist of different materials or are doped differently. However, when the device calculation is set up, the right electrode Fermi level is shifted (and all energy eigenvalues along with it) such that the Fermi levels coincide at zero bias.

The current is independent of how the voltages are applied, it only depends on their difference (i.e. on \(V_{bias}\)).

A positive current flows when the applied voltage for the left electrode is higher than that of the right one, thus the positive current direction is left to right.

The TransmissionSpectrum itself is usually rather insensitive to the temperatures \(T_L\) and \(T_R\) used in the self-consistent calculation, but often depends strongly on the electrode voltages \(V_L\) and \(V_R\). Thus, for an accurate estimate of the current, the TransmissionSpectrum should be calculated self-consistently for each desired bias. It is, however, possible to obtain an approximation for the current at a different bias than that used in the self-consistent calculation by giving a new set of voltages to the current() method. This will simply modify the integration range on the existing transmission spectrum. The temperature-dependence of the current can, on the other hand, in most cases be rather accurately estimated by only changing the electrode temperatures in the current calculation.

Differential conductance

The differential conductance is calculated from the transmission spectrum using:

\[\sigma(V_L, V_R, T_L, T_R, \alpha_L, \alpha_R) =\lim_{\delta V \rightarrow 0} \frac{I(V_L+\alpha_L \delta V, V_R- \alpha_R \delta V,T_L, T_R )}{\delta V}.\]

The coupling constants \(\alpha_L + \alpha_R = 1.0\) models how the transmission spectrum couples with the left and right electrode. For a molecule strongly bound to the left electrode, we must have \(\alpha_L=0\) and \(\alpha_R=1\).

A more accurate estimate of the differential conductance which avoid introducing \(\alpha_L\) and \(\alpha_R\) is obtained by calculating the self-consistent current at a number of applied biases, \(V_{bias}^1, V_{bias}^2, \ldots\), and performing numerical differentiation:

\[\sigma(V_{bias}, T_L, T_R) = \frac{I(V_{bias}^1,T_L, T_R)-I(V_{bias}^2,T_L, T_R))}{V_{bias}^1-V_{bias}^2}.\]

Default automatic energy window

If no energies are provided, TransmissionSpectrum attempts to estimate an energy range such that the whole bias window is covered with an accurate enough spacing between the energy points for the integration.

The lower and upper bounds for the automatic energy window are defined as follows. First, we define an energy padding on each side of the bias window as \(E^p_L = 30 k_B T_L\) and \(E^p_R = 30 k_B T_R\), for the left and right electrodes, respectively. This is to ensure that the tails of the Fermi distributions are included in the energy window with a very tight accuracy. Then, we use the following formula to calculate the minimum and maximum energies of the energy window:

\[ \begin{align}\begin{aligned}E_{\mathrm{lower}} = \min(E_{F, L} - E^p_L, E_{F, R} - E^p_R) - E_0,\\E_{\mathrm{upper}} = \max(E_{F, L} + E^p_L, E_{F, R} + E^p_R) - E_0,\end{aligned}\end{align} \]

where \(E_0\) denotes the energy zero as specified by the energy_zero_parameter, and \(E_{F, L}\) and \(E_{F, R}\) denote the Fermi levels of the left and right electrodes, respectively. The spacing between the energy points in the energy window is determined by the real-axis point density of the non-equilibrium contour (\(\delta E_{\mathrm{neq}}\)) used for the calculation (see Complex contour integration for more details about the contour integration). However, we always ensure that there are at least three points per each \(k_B T\) value for the lower electrode temperature, i.e. the energy spacing can be written as

\[\delta E = \min(\delta E_{\mathrm{neq}}, \frac{\min(k_B T_L, k_B T_R)}{3}).\]

As an example, for electrode temperatures \(T_L =\) 300 K, \(T_R =\) 400 K and a default RealAxisContour as the non-equilibrium contour, we obtain the (absolute) energy window [-1.084 eV, 0.984 eV] with 241 equally spaced energy points.