STM images

The STM is a revolutionary experimental surface probe that has provided direct local insight into the surface electronic structure. Sometimes the interpretation of STM topographs are not straightforward and therefore theoretically modeled STM images may resolve conflicting possibilities and point to an underlying atomistic model. ASE includes python modules for generating Tersoff-Hamann STM topographs.

The calculated tunneling current will be proportional to:

\[\int_{\epsilon_F}^{\epsilon_F+eV} \sum_{kn} w_{\mathbf k} |\Psi_{\mathbf k n}(\mathbf r)|^2 \delta(\epsilon - \epsilon_{\mathbf k n}) d\epsilon,\]

where \(V\) is the bias voltage, \(w_{\mathbf k}\) is the \(\mathbf k\)-point weight and \(\Psi_{\mathbf k n}(\mathbf r)\) is the wave function.

More details:

class ase.dft.stm.STM(atoms, symmetries=None, use_density=False)[source]

Scanning tunneling microscope.

atoms: Atoms object or filename

Atoms to scan or name of file to read LDOS from.

symmetries: list of int

List of integers 0, 1, and/or 2 indicating which surface symmetries have been used to reduce the number of k-points for the DFT calculation. The three integers correspond to the following three symmetry operations:

[-1  0]   [ 1  0]   [ 0  1]
[ 0  1]   [ 0 -1]   [ 1  0]
use_density: bool

Use the electron density instead of the LDOS.

calculate_ldos(bias)[source]

Calculate local density of states for given bias.

find_current(ldos, z)[source]

Finds current for given LDOS at height z.

get_averaged_current(bias, z)[source]

Calculate avarage current at height z (in Angstrom).

Use this to get an idea of what current to use when scanning.

linescan(bias, current, p1, p2, npoints=50, z0=None)[source]

Constant current line scan.

Example:

stm = STM(...)
z = ...  # tip position
c = stm.get_averaged_current(-1.0, z)
stm.linescan(-1.0, c, (1.2, 0.0), (1.2, 3.0))
pointcurrent(bias, x, y, z)[source]

Current for a single x, y, z position for a given bias.

scan(bias, current, z0=None, repeat=(1, 1))[source]

Constant current 2-d scan.

Returns three 2-d arrays (x, y, z) containing x-coordinates, y-coordinates and heights. These three arrays can be passed to matplotlibs contourf() function like this:

>>> import matplotlib.pyplot as plt
>>> plt.contourf(x, y, z)
>>> plt.show()
scan2(bias, z, repeat=(1, 1))[source]

Constant height 2-d scan.

Returns three 2-d arrays (x, y, I) containing x-coordinates, y-coordinates and currents. These three arrays can be passed to matplotlibs contourf() function like this:

>>> import matplotlib.pyplot as plt
>>> plt.contourf(x, y, I)
>>> plt.show()
sts(x, y, z, bias0, bias1, biasstep)[source]

Returns the dI/dV curve for position x, y at height z (in Angstrom), for bias from bias0 to bias1 with step biasstep.

write(filename)[source]

Write local density of states to JSON file.