[1]:
from ase.io import read
from abtem.waves import PlaneWave
from abtem.waves.transfer import CTF
import matplotlib.pyplot as plt
WARNING:hyperspy_gui_traitsui:The module://matplotlib_inline.backend_inline matplotlib backend is not compatible with the traitsui GUI elements. For more information, read http://hyperspy.readthedocs.io/en/stable/user_guide/getting_started.html#possible-warnings-when-importing-hyperspy.
WARNING:hyperspy_gui_traitsui:The traitsui GUI elements are not available.
Contrast Transfer Function¶
The Contrast Transfer Function (CTF) describes the aberrations of the objective lens in HRTEM and specifies 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.
Polar expansion of the phase error¶
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. This deviation can be represented as a phase error, \(\chi(k, \phi)\), of the spatial frequency in polar coordinates. The phase error can be written as a series expansion
where \(n,m\) index the different aberrations.
If the microscope is well aligned then off-axis aberrations (astigmatisms) are small and the phase error is dominated by the first two isotropic terms
where \(\Delta f = -C_{1,1}\) is the defocus and \(C_s=C_{3,3}\) is the third order spherical aberration.
The CTF
object takes parameters of the form Cnm
and phinm
, and may also be given using their common aliases, e.g. defocus = -C10
. The expansion is implemented up to 5th order, and all the coefficients labels are shown in the table below.
|
|
|
||||
|
|
|
|
|||
|
|
|
|
|
||
|
|
|
|
|
|
|
|
|
|
|
|
|
|
abTEM exclusively uses polar expansion coefficients; however, a conversion utility from the Cartesian to the polar representation is available in abtem.transfer.cartesian2polar
.
Creating a Contrast Transfer Function and applying phase aberrations¶
We create a CTF
object with an electron energy of \(300 \ \mathrm{keV}\), a spherical aberration of \(-7~\mu \mathrm{m}\) (remember that abTEM uses units of Å) and the Scherzer defocus (here \(-45.5 \ \mathrm{Å}\)).
[2]:
Cs = -7e-6 * 1e10
ctf = CTF(Cs=Cs, energy=300e3)
ctf.defocus = ctf.scherzer_defocus
print(f"CTF defocus: {ctf.defocus:.2f} mrad")
CTF defocus: -45.47 mrad
This can also be accomplished by setting the parameters with a dictionary.
[3]:
aberrations = {'C10': -ctf.defocus, 'C30': Cs}
ctf2 = CTF(aberrations=aberrations, energy=300e3)
We can use the method .show
to preview the CTF along a radial direction up to a given angle (in milliradians).
[4]:
ctf2.show(angular_units=True)
[4]:
<AxesSubplot:xlabel='x [1 / Å]'>

[ ]:
#ctf2.show(max_semiangle=40);
Given the exit wave function \(\phi_{exit}\) and the phase error \(\chi\), the wave function at the image plane is given by
To perform this operation, we use the apply_ctf
method. Here we apply the CTF to the exit wave from an earlier tutorial.
[5]:
exit_wave = PlaneWave(energy=300e3, sampling=.05).multislice(read('data/srtio3_110.cif'))
image_wave1 = exit_wave.apply_ctf(ctf)
image_wave1.show();

Aperture and partial coherence¶
The CTF is further limited by any aperture present and by partial coherence. This is conventionally described as a multiplication with the aperture function, \(A(k)\), and a partial coherence envelope, \(E(k)\). Hence, the wave function at the image plane becomes
The aperture function cuts off beams scattered above a certain critical angle
We will cut off the CTF at the angle corresponding to the Scherzer point resolution, which is defined as the angle (in radians) where the CTF crosses the abscissa for the first time (here converted to mrad by the factor 1000).
[6]:
ctf.semiangle_cutoff = ctf.crossover_angle
print(f"Semi-angle cutoff: {ctf.crossover_angle:.2f} mrad")
image_wave = exit_wave.apply_ctf(ctf)
fig,(ax1,ax2) = plt.subplots(1,2, figsize=(9,4))
ctf.show(50, ax=ax1)
image_wave.show(ax=ax2);
Semi-angle cutoff: 36.04 mrad

[8]:
ctf.fourier_space_image(exit_wave).show(max_angle=70, angular_units=True)
[8]:
(<AxesSubplot:xlabel='Scattering angle x [mrad]', ylabel='Scattering angle y [mrad]'>,
<matplotlib.image.AxesImage at 0x1c5887f2e20>)

The partial coherence envelope function \(E(k)\) smoothly dampens the signal of beams scattered to high angles. Its most important contribution is generally the focal spread, or, equivalently the chromatic aberration.
[9]:
ctf.focal_spread = 30
image_wave2 = exit_wave.apply_ctf(ctf)
fig,(ax1,ax2) = plt.subplots(1, 2, figsize=(9,4))
ctf.show(50, ax=ax1)
image_wave2.show(ax=ax2);

Probes and CTFs¶
The contrast transfer function can also be calculated for a focused probe.
[26]:
from abtem.waves import Probe
probe = Probe(sampling=.025, extent=10, energy=80e3, semiangle_cutoff=30, defocus=100)
probe.ctf.show()
[26]:
<AxesSubplot:xlabel='x [1 / Å]'>

We can further model the effect of different types of aberrations on our probe by changing the parameters of the contrast transfer function.
[27]:
fig, (ax1, ax2, ax3) = plt.subplots(1, 3, figsize=(12,4))
probe_copy = probe.copy()
probe_copy.aberrations.set_parameters({'C12': 20, 'phi12': 0.785})
probe_copy.show(ax=ax1)
ax1.set_title('astigmatism')
probe_copy = probe.copy()
probe_copy.aberrations.set_parameters({'C21': 5000, 'phi21': 0.785})
probe_copy.show(ax=ax2)
ax2.set_title('coma')
probe_copy = probe.copy()
probe_copy.aberrations.set_parameters({'C23': 2500, 'phi23': 0.785})
probe_copy.complex_images().show(ax=ax3);
ax3.set_title('three-fold astigmatism')
[27]:
Text(0.5, 1.0, 'three-fold astigmatism')

Partial coherence¶
The aperture function cuts off beams scattered above a certain critical angle, however it is often the partial coherence that effectively dampens the signal of beams scattered at high angles, and imposes a maximum to the transmitted spatial frequency.
In abTEM the envelope function is constructed as a product of three envelopes
Below we create a CTF
object with all the sources of partial coeherence.
[9]:
from abtem import CTF
ctf = CTF(energy=200e3, defocus=600, Cs=1e7, focal_spread=100, angular_spread=0.5, gaussian_spread=2)
ctf.show(30);

Note: The quasi-coherent approximation
Partial spatial coherence is in the contrast transfer function using the quasi-coherent approximation, where high spatial frequencies are damped using an envelope function. This approximation assumes that angular and focal spread may be summed coherently. While this is efficient, the approximation may not always be appropriate.
Partial temporal coherence¶
A small spread in energy, \(\Delta E\), of the incident electrons is due to the chromatic aberration of the objective lens equivalent to a small spread in defocus. Fluctuations in the focusing currents, \(\Delta I\), of the objective lens also produce an incoherent spread in defocus. Combining these effects, the \(1 / e\) width (here \(e \approx 2.71828\) denotes the Euler number) of the distribution of defocus values (the focal spread) can be written as
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, it can be shown that the focal spread can be approximated as
The parameter \(\delta\) is equivalent to focal_spread
.
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 the angular spread. Assuming that each incident direction performs its own scattering and that the distribution of incident directions is small, then it can be shown that the angular spread can be modelled by the spatial coherence envelope function, as given by
where \(\beta\) is the \(1/e\) width of the distribution of angles. The parameter \(\beta\) is equivalent to angular_spread
.
Deflections and thermal magnetic field noise¶
Blurring can also be caused by all kinds of noise leading to a random deflection of the image relative to the detector, such as vibrations, drift of the stage, and magnetic noise fields resulting from eddy currents in the material of the lenses. Assuming that the image deflection is small and follows a Gaussian distribution, this can be included as an additional envelope
where \(\sigma\) is the \(1/e\) width of the deflection. The parameter \(\sigma\) is equivalent to the gaussian_spread
.