Hide code cell source
import ase
import matplotlib.pyplot as plt
import numpy as np
from abtem.transfer import polar_symbols

import abtem

abtem.config.set({"local_diagnostics.progress_bar": False})
<abtem.core.config.set at 0x1e3ea507d10>

Contrast transfer function#

The contrast transfer function (CTF) describes the aberrations of the objective lens in HRTEM and how the condenser system shapes the probe in STEM. Here we describe how to create a CTF with specific aberrations and how this affects the resulting images.

Phase aberrations#

Ideally, a lens forms a spherical wave converging on or emerging from a single point. In practice, aberrations cause the wave front to deviate from a spherical surface. In a polar representation, this deviation, or phase aberration, can be written as \(\chi(k, \phi)\), where \(k = \sqrt{k_x^2 + k_y^2}\) is the radial component and \(\phi = \arctan(k_y / k_x)\) is the azimuthal component. We expand the phase error as a series expansion

\[ \chi(k, \phi) = \frac{2 \pi}{\lambda} \sum_{n,m} \frac{1}{n + 1} C_{n,m} (k \lambda)^{n+1} \cos\left[m (\phi - \phi_{n,m}) \right] \quad , \]

where \(n\) and \(m\) are the radial and azimuthal order.

If the microscope is well aligned then off-axis aberrations (astigmatisms and comas) are small and the aberrations are dominated by the first two isotropic terms. The aberrations reduce to an expression you may be more familiar with

\[ \chi(k) \approx \frac{2\pi}{\lambda}\left( \frac{\lambda^2 k^2}{2} \Delta f + \frac{\lambda^4 k^4}{4} C_s \right) \quad , \]

where \(\Delta f = -C_{1,0}\) is the defocus and \(C_s=C_{3,0}\) is the third order spherical aberration.

In abTEM, the CTF object takes parameters of the form Cnm and phinm, and some may alternatively also be given using their common aliases, e.g. defocus = -C10. The expansion is implemented up to 5th order and all the coefficients keywords and their aliases are given in the tables below.

\(m = 0\)

\(m = 1\)

\(m = 2\)

\(m = 3\)

\(m = 4\)

\(m = 5\)

\(m = 6\)

\(n = 1\)

C10

C12, phi12

\(n = 2\)

C21, phi21

C23, phi23

\(n = 3\)

C30

C32, phi32

C34, phi34

\(n = 4\)

C41, phi41

C43, phi43

C45, phi45

\(n = 5\)

C50

C52, phi52

C54, phi54

C56, phi56

\(m = 0\)

\(m = 1\)

\(m = 2\)

\(m = 3\)

\(m = 4\)

\(m = 5\)

\(m = 6\)

\(n = 1\)

defocus

astigmatism

\(n = 2\)

coma

trefoil

\(n = 3\)

Cs

astigmatism2

quadrafoil

\(n = 4\)

coma2

trefoil2

pentafoil

\(n = 5\)

C5

astigmatism3

quadrafoil2

hexafoil

Below we visualize each aberration by mapping the phase to different colors. We can also visualize each phase aberration by showing its effect on a small probe wave function. This is called the point spread function (PSF) and describes the response of an imaging system to a point source. Note that the magnitude of the aberrations are scaled by a power law are scaled by the radial order (\(n\)) according to a power law.

ctf

ctf

Applying phase aberrations to wave functions#

Given an exit wave function \(\phi_{exit}\) and an objective lens with a phase error \(\chi\) (and infinite aperture), the wave function at the image plane is given by

(4)#\[ \psi_{image}(k, \phi) = \psi_{exit}(k, \phi) \exp(-i \chi(k, \phi)) \quad . \]

We start by simulating an exit wave for an MoS2 model. An atomic model with a hexagonal cell is created using the ASE’s mx2 method, which we make orthogonal using orthogonalize_cell (see our walkthrough on atomic models).

atoms = ase.build.mx2(vacuum=2)

atoms = abtem.orthogonalize_cell(atoms)

fig, (ax1, ax2) = plt.subplots(1, 2, figsize=(7, 4))
abtem.show_atoms(atoms, ax=ax1, title="Beam view", legend=True)
abtem.show_atoms(atoms, ax=ax2, plane="xz", title="Side view");
../../_images/c9f763e1d4add055e611bf2b7ba2b5128b80c84bb2a2b653a280938bb889e8c4.png

We run the multislice algorithm for an incident plane wave at \(80 \ \mathrm{keV}\) through a default potential (using the Lobato indepdent atom parameterization, see our walkthrough on potentials).

exit_wave = abtem.PlaneWave(energy=80e3, sampling=0.02).multislice(atoms)

We create a compatible CTF for electrons of the same energy of \(80 \ \mathrm{keV}\). The spherical aberration will be set to \(-20~\mu \mathrm{m}\) (remember that abTEM uses units of \(\mathrm{Å}\)) and the defocus is set to the Scherzer defocus (here \(-111.92 \ \mathrm{Å}\)).

Cs = -20e-6 * 1e10  # 20 micrometers

ctf = abtem.CTF(Cs=Cs, energy=80e3, defocus="scherzer")

print(f"defocus = {ctf.defocus:.2f} Å")
defocus = -111.92 Å

The aberrations may also be given as a dictionary.

aberration_coefficients = {"C10": -ctf.defocus, "C30": Cs}

ctf = abtem.CTF(aberration_coefficients=aberration_coefficients, energy=ctf.energy)

We can get radial LineProfiles up to a given angle to preview the CTF, which by convention display the imaginary part of complex exponential of the phase error.

visualization = ctf.profiles(max_angle=60).show(units="mrad")
../../_images/ca82bc1913b4db7a0bd33532a7e144fd7c606e538fb3abbed851e773f63e8c15.png

To perform the multiplication in Eq. (4), we use the apply_ctf method.

image_wave = exit_wave.apply_ctf(ctf)

We show the exit wave function and the wave function at image plane below.

abtem.stack([exit_wave, image_wave], ("Exit wave", "Aberrated wave")).show(
    explode=True,
    figsize=(14, 5),
    common_color_scale=True,
    cbar=True,
);
../../_images/d4b1a0d6793ba855a78376c67f6bd342d3ea76d8fdd8aada3bc4d9912fbd5ac2.png

Without an aperture, the high-frequency components of the CTF lead to the observed complicated contrast variation.

Focal series#

It common that the aberrations of your microscope are not known exactly, and thus it is useful to simulate an ensemble of images with different aberrations.

We create an ensemble of CTFs for an range of defocii.

ctf_ensemble = abtem.CTF(Cs=-25e-6 * 1e10, energy=80e3)

ctf_ensemble.C10 = np.linspace(0, -ctf_ensemble.scherzer_defocus, 5)

When this CTF is applied to the wave function we obtain an ensemble of wave functions.

waves_image_series = exit_wave.apply_ctf(ctf_ensemble)

waves_image_series.show(
    figsize=(14, 5),
    cbar=True,
    explode=True,
    common_color_scale=True,
);
../../_images/2630f0b4c34762655305048d0c0b885e49a143f60b9a5bea80315b24649b5b09.png

Aperture#

In practice, the maximum transferred frequency is limited by the aperture of the objective lens. This is conventionally described as a multiplication with the aperture function:

\[ \psi_{\mathrm{image}}(k, \phi) = A(k) \exp[-i \chi(k, \phi)] \psi_{\mathrm{exit}}(k, \phi), \]

where \(A(k)\) is the aperture function

\[\begin{split} A(k) = \begin{cases} 1 & \text{if } x \leq k_{cut} \\ 0 & \text{if } x > k_{cut} \end{cases} . \end{split}\]

We will cut off the CTF at the angle corresponding to the Scherzer point resolution, which is defined as the angle where the phase of the CTF crosses the abscissa for the first time (crossover_angle).

ctf_cutoff = ctf.copy()
ctf_cutoff.semiangle_cutoff = ctf_cutoff.crossover_angle

print(f"cutoff = {ctf_cutoff.crossover_angle:.3f} mrad")
cutoff = 33.455 mrad
ctf_cutoff.profiles(max_angle=40).show(legend=True, units="mrad");
../../_images/c34455b9a373ffcb542bdeb87286545141f6df8715356cb636fd664e423b3a42.png

We show the resulting image below and see, as expected, that applying this aperture removes high-frequency information from the wave function intensity.

waves_cutoff = exit_wave.apply_ctf(ctf_cutoff)

abtem.stack(
    [exit_wave, image_wave, waves_cutoff], ("Exit wave", "Imaged wave", "+ cutoff")
).show(
    explode=True,
    figsize=(14, 5),
    common_color_scale=True,
    cbar=True,
);
../../_images/5dff945eec933504ab8066b347d1458cecbb3151eed629986ed85e633615f735.png

Asymmetric aberrations#

Assymetric aberrations such as astigmatism also requires setting an angle, here we adjust the CTF by adding \(50 \ \mathrm{Å}\) astigmatism (\(C_{12}\)) at angles from \(0\) to \(\pi / 4 \ \mathrm{rad}\).

Note

abTEM exclusively uses polar expansion coefficients; however, a conversion utility from the Cartesian to the polar representation is available in abtem.waves.transfer.cartesian2polar. Note that there several are different conventions for aberration coefficients, abTEM uses \(\mathrm{\AA}\) for all amplitude parameters and radians for the azimuthal parameters.

ctf_asymmetric = ctf_cutoff.copy()

ctf_asymmetric.astigmatism = 50
ctf_asymmetric.astigmatism_angle = np.linspace(0, np.pi / 2, 3)

We can no longer simply visualize the aberrations as a line profile, as the profile will depend on the azimuthal angle. We can instead evaluate the CTF as an ensemble of complex reciprocal space images and visualize those.

measurements = ctf_asymmetric.to_diffraction_patterns(
    max_angle=ctf_cutoff.crossover_angle * 1.01
)
vis = measurements.show(
    explode=True, figsize=(12, 4), cbar=True, common_color_scale=True
)
../../_images/d3a6b38434592e91f2b6222d8358ffbda28e04e66cabd6667e82a6ad183a3f6b.png

Alternatively, we can get the point spread functions as complex real-space images.

psf = ctf_asymmetric.to_point_spread_functions(gpts=128, extent=15).compute()

psf.show(explode=True, figsize=(12, 4), cbar=True, common_color_scale=True, vmax=0.51e-3, cmap="hsluv");
../../_images/0ce67b1a99fb314089612cda15a99a999df0a869b7de7426f0fa31811d7204d8.png

Partial coherence (quasi-coherent)#

Partial coherence acts similarly to the aperture function to dampen the high spatial frequencies of the signal. Partial coherence may be approximated by multiplying the CTF by envelope functions:

\[ \psi_{\mathrm{image}}(k, \phi) = \mathrm{CTF}(k, \phi) \psi_{\mathrm{exit}}(k, \phi), \]

where the \(\mathrm{CTF}\) is now given as

\[ \mathrm{CTF}(k, \phi) = E_t(k) E_s(k) A(k) \exp[-i \chi(k, \phi)] \]

and \(E_t(k)\) and \(E_s(k)\), to be defined below, are the temporal and spatial envelopes, respectively.

Warning

Here partial coherence is applied using the convenient quasi-coherent approximation. This is only appropriate for HRTEM simulations when the effect is sufficiently small, but it is generally never appropriate for STEM simulations. The more accurate approach of using an incoherent summation is introduced in a tutorial.

Partial temporal coherence#

The most important source of partial coherence in HRTEM is partial temporal coherence. A small spread in energy of the incident electrons, \(\Delta E\), is due to the chromatic aberration of the objective lens equivalent to a small spread in defocus. Fluctuations in the focusing currents of the objective lens, \(\Delta I\), also produce an incoherent spread in defocus. Combining these effects, the standard deviation of the focal spread may be written as

\[ \delta = C_c \sqrt{4 \left(\frac{\Delta I_\text{obj}}{I_\text{obj}}\right)^2 + \left(\frac{\Delta E}{V_\text{acc}}\right)^2 + \left(\frac{\Delta V_\text{acc}}{V_\text{acc}}\right)^2} \quad. \]

The terms \(\Delta I_\text{obj}/I_\text{obj}\) and \(\Delta V_\text{acc}/V_\text{acc}\) represent instabilities in of the total current in the magnetic lenses and the acceleration voltage. \(\Delta E/V_\text{acc}\) is the energy spread of electrons emitted by the source. \(C_c\) is the chromatic aberration coefficient. Assuming \(\delta\) is small, and focal spread follows a Gaussian distribution, the temporal envelope may be written [Kir10]

\[ E_t(k) = \exp\left[-\frac{1}{4} \left(\pi \lambda \delta \right)^2 k^4 \right] \quad . \]

The parameter \(\delta\) above is equivalent to the focal_spread keyword and given in \(\mathrm{Å}\).

We calculate a realistic value of the focal spread to be \(52.5 \ \mathrm{Å}\). This was estimated using a typical energy spread of \(0.35 \ \mathrm{eV}\) (standard deviation) for a non-monochromated beam and a typical chromatic aberration of \(1.2 \ \mathrm{mm}\) which would be considered good for a microscope without chromatic aberration correction. We assume that the contributions to the focal spread from instabilities in the magnetic lenses are negligible.

Cc = 1.2e-3 * 1e10
energy = exit_wave.energy
energy_spread = 0.35

focal_spread = Cc * energy_spread / energy

print(f"focal spread = {focal_spread:.2f} Å")
focal spread = 52.50 Å

We show a plot of the resulting CTF below. Given this we may conclude that for a focal spread of this magnitude, the aperture will have almost no effect as the envelope already goes to zero at the cutoff.

ctf_partial = ctf_cutoff.copy()
ctf_partial.focal_spread = focal_spread
ctf_partial.profiles().show(units="mrad", legend=True);
../../_images/e7d37126b25955aab3db2cf3a7bcdb799bf6bbd4720a7935c954e0ccf33c5fa6.png

Partial spatial coherence#

As the electron source has a finite size, the incident beam contains a distribution of incident directions. In HRTEM this is quantified by angular spread. Assuming that each incident direction performs its own scattering and that the distribution of incident directions is small, it can be shown that the angular spread can be modelled by the spatial coherence envelope function[Kir10], as given by

\[ E_s(k) = \exp\left(-\frac{\beta}{4\lambda^2}\left| \frac{\partial \chi(k)}{\partial k}\right|^2 \right) \quad . \]

where \(\beta\) is the \(1/e\) width of the distribution of angles. The parameter \(\beta\) is equivalent the angular_spread keyword and given in \(\mathrm{mrad}\).

Given a modern electron source such as a field emission gun, the partial spatial coherence is typically negligible, even in microscopes without chromatic aberration correction.

Given an angular spread of \(1 \ \mathrm{mrad}\), which would be considered large, we obtain the spatial envelope illustrated below, which we see has almost no effect in the presence of partial temporal coherence.

ctf_partial.angular_spread = 1

ctf_partial.profiles().show(legend=True);
../../_images/8ffd96985bb12a8f0c6fb9681f4c1f92a0baef4c4486bd2f5bb37c3a064766b1.png

Finally we apply the CTF with the partial coherence envelope to the wave function. The effect is to further blur the image and to lower the contrast.

wave_partial = exit_wave.apply_ctf(ctf_partial)

abtem.stack(
    [exit_wave, image_wave, waves_cutoff, wave_partial],
    ("Exit wave", "Imaged wave", "+ cutoff", "+ partial coherence"),
).show(
    explode=True,
    figsize=(14, 5),
    common_color_scale=True,
    cbar=True,
);
../../_images/15010cc54ce3f93456a68f35b90a4882333025cfcf89dee36e088e5029356a76.png

Note

Blurring can also be caused by noise leading to random deflection of the image relative to the detector, such as vibrations, drift of the stage and magnetic noise in the lenses. Such effects may be reasonably included by applying a Gaussian blur.