[1]:
import numpy as np
from ase.io import read
from abtem import *

Scan and detect

Scanning imaging modes such as STEM works by rastering an electron probe across a sample pixel by pixel, and recording the scattering signal. As in earlier tutorials, we start by importing our graphene atoms and creating a probe, here at 80 keV.

[2]:
atoms = read('data/orthogonal_graphene.cif')

probe = Probe(energy=80e3, semiangle_cutoff=30, focal_spread=60, defocus=50)

potential = Potential(atoms, sampling=.04)

Line scan

The LineScan function scans the probe along a straight line between two points with a spacing given by the gpts parameter.

Below we create a linescan object to scan between two atoms, and run a multislice simulation to get the corresponding exit wave. The scan objects can be given as an argument of the plotting method of 2D arrays or the show_atoms function to visualise the scan region.

[3]:
from abtem.scan import LineScan

linescan = LineScan(start=[2 * np.sqrt(3) * 1.42, 0], end=[2 * np.sqrt(3) * 1.42, 3 * 1.42], gpts=40, endpoint=False)

ax, im = potential.project().show()
linescan.add_to_mpl_plot(ax)
../_images/walkthrough_06_scan_and_detect_6_0.png
[4]:
line_exit_probes = probe.multislice(positions=linescan.get_positions(), potential=potential)

Detect

In real experiments, the exit probes are not directly detected. Correspondingly, in abTEM the way a signal is detected is controlled by detector objects. There are several different types of detectors, but for now we focus on the basic AnnularDetector, which depending on the choice of integration region can represent the detector used in bright-field, medium- or high-angle annular dark-field microscopy. The integration region is given by an inner and outer radius in mrad; below we create three different types of detectors, and plot their respective integration regions.

[5]:
from abtem.detect import AnnularDetector
import matplotlib.pyplot as plt

bright = AnnularDetector(inner=0, outer=20)
haadf = AnnularDetector(inner=90, outer=200)

fig, (ax1, ax2) = plt.subplots(1, 2, figsize=(8,4))

bright.show(probe, ax=ax1, title='Bright field')
haadf.show(probe, ax=ax2, title='HAADF');
../_images/walkthrough_06_scan_and_detect_10_0.png

The probes can be detected using the .detect method, this condenses each probe down to a single integrated intensity.

[6]:
haadf.detect(line_exit_probes)
[6]:
array([2.01043309e-04, 1.96584660e-04, 1.84240096e-04, 1.66724407e-04,
       1.47861909e-04, 1.31806446e-04, 1.22136480e-04, 1.21057165e-04,
       1.28894375e-04, 1.44012854e-04, 1.63181350e-04, 1.82292832e-04,
       1.97260582e-04, 2.04878073e-04, 2.03432137e-04, 1.92936321e-04,
       1.74967878e-04, 1.52167384e-04, 1.27575055e-04, 1.03979575e-04,
       8.34414677e-05, 6.70754162e-05, 5.50919103e-05, 4.70270970e-05,
       4.20563701e-05, 3.92923321e-05, 3.80103447e-05, 3.77825891e-05,
       3.85398089e-05, 4.05770988e-05, 4.45109981e-05, 5.11726430e-05,
       6.14163800e-05, 7.58523311e-05, 9.45333741e-05, 1.16693598e-04,
       1.40634496e-04, 1.63840887e-04, 1.83366632e-04, 1.96421548e-04],
      dtype=float32)

The typical syntax for simulating a STEM experiment requires giving the detector(s) as an argument of the scan object.

[7]:
measurements = probe.scan(linescan, [haadf, bright], potential, pbar=True)


The output is a list of measurements with one entry per detector, e.g. the bright-field one.

[8]:
measurements[0]
[8]:
<abtem.measure.Measurement at 0x19618170f40>
[9]:
fig, (ax1, ax2) = plt.subplots(1,2, figsize=(12, 4))
measurements[0].show(ax=ax1, title='Bright field')
measurements[1].show(ax=ax2, title='HAADF');
../_images/walkthrough_06_scan_and_detect_17_0.png

We can write a measurement to disk into an .hdf5 file.

[10]:
measurements[0].write('data/linescan.hdf5')
[10]:
'data/linescan.hdf5'

Grid scan

The GridScan object is used essentially the same as the LineScan object, except the start and end coordinates now denote the corners of the grid. In the following, we scan over a tileable section of the potential. By specifying an output file for the detector, we write the measurements directly to disk as they are made.

We set the probe step size (or sampling) slightly better than the Nyquist frequency of the probe contrast transfer function. The resulting image will be quite pixelated, however, we can use interpolation to improve the sampling.

[11]:
from abtem.scan import GridScan

detector = AnnularDetector(inner=86, outer=190, save_file='data/gridscan.hdf5')
gridscan = GridScan(start=[0, 0], end=[np.sqrt(3) * 1.42, 3 * 1.42], sampling=probe.ctf.nyquist_sampling * .9)

ax, im = potential.project().show();

gridscan.add_to_mpl_plot(ax)
../_images/walkthrough_06_scan_and_detect_21_0.png

The number of positions in a densely sampled grid can obviously be large and the calculation take some time, hence it may be useful to show the progress as below.

[12]:
measurement_files = probe.scan(gridscan, [detector], potential, pbar=True)


We read back the measurements from disk and tile the image.

[13]:
from abtem.measure import Measurement

measurement = Measurement.read('data/gridscan.hdf5')
measurement.show();
../_images/walkthrough_06_scan_and_detect_25_0.png

The pixel size in the image is set by our sampling, which affects the computational cost of the calculation. The resulting images can be tiled and smoothed by interpolating values between the pixels (here with a virtual sampling of \(0.04\) Å).

[14]:
new_measurement = measurement.tile((7, 4))
new_measurement = new_measurement.interpolate(.04)

new_measurement.show();
../_images/walkthrough_06_scan_and_detect_27_0.png

It is not possible to save to disk during the scan, but it is straightforward to export the measurement as an image.

[15]:
new_measurement.save_as_image('data/gridscan.tif')