Hide code cell source
%config InlineBackend.rc = {"figure.dpi": 72, "figure.figsize": (6.0, 4.0)}
%matplotlib inline

import ase
import matplotlib.pyplot as plt
from ase.io import read

import abtem

SAED quickstart#

This notebook demonstrates a basic simulation of selected area electron diffraction of SiO2 zeolite.

Configuration#

We start by (optionally) setting our configuration. See documentation for details.

abtem.config.set({"device": "cpu", "fft": "fftw"})
<abtem.core.config.set at 0x1ccfe5cfd90>

Atomic model#

We import a unit cell from a .cif-file. See our walkthough or our tutorial on atomic models.

sio2 = read("data/SiO2.cif")

abtem.show_atoms(sio2, plane="xy", legend=True);
../../../_images/01599072b2cd671ab3904f382eba8d027749685a662f89e44aee4a0f6202c6f4.png

The structure is rotated such that the previous \(xy\)-plane is rotated into the \(xz\)-plane. The structure is repeated along \(z\) to obtain the required thickness. If frozen phonons are required the structure should also be repeated in the \(xy\) to describe a thermal ensemble.

sio2_rotated = abtem.atoms.rotate_atoms_to_plane(sio2, "xz")

sio2_repeated = sio2_rotated * (3, 2, 20)

fig, (ax1, ax2) = plt.subplots(1, 2, figsize=(12, 4))
abtem.show_atoms(sio2_repeated, ax=ax1, plane="xy", title="Beam view")
abtem.show_atoms(sio2_repeated, ax=ax2, plane="yz", title="Side view", linewidth=0);
../../../_images/05c06fdafde66aec77988d71bcfa8fa0854983526cf039c9330d1efeb49c5d91.png

Potential#

We create an ensemble of potentials using the frozen phonon model. See our walkthrough on frozen phonons.

frozen_phonons = abtem.FrozenPhonons(sio2_repeated, 8, sigmas=0.1)

We create a potential from the frozen phonons model, see walkthrough on potentials.

potential = abtem.Potential(
    frozen_phonons,
    gpts=512,
    projection="infinite",
    slice_thickness=2,
    exit_planes=50,
)

Wave function#

We create a plane wave function at an energy of \(200 \ \mathrm{keV}\). See our walkthrough on wave functions.

wave = abtem.PlaneWave(energy=200e3)

Multislice#

We run the multislice algorithm to calculate the exit waves, see our walkthrough on multislice.

exit_waves = wave.multislice(potential)

We calculate the diffraction patterns up to a scattering angle of \(10 \ \mathrm{mrad}\).

measurement_ensemble = exit_waves.diffraction_patterns(max_angle=10)

measurement_ensemble.shape
(8, 5, 35, 35)

The result is an ensemble of images, one for each frozen phonon and exit plane, we average the ensemble across the phonon dimensions obtain the final thickness series.

measurement = measurement_ensemble.mean(0)

measurement.compute()
[########################################] | 100% Completed | 15.56 s
<abtem.measurements.DiffractionPatterns object at 0x000001CC8C277370>

Visualize results#

We show the thickness series as an exploded plot. Before plotting the direct beam is blocked as it typically has a much higher intensity than the scattered beams.

visualization = measurement.block_direct().show(
    explode=True,
    figsize=(18, 5),
    cbar=True,
    common_color_scale=True,
)
../../../_images/f044953987d4d5743ae1a7860555427bd650c18834a7a0989524ed302dbd2bc0.png

Index diffraction patterns#

Instead of the pixelated representation above, the diffraction patterns can be represented as a collection indexed diffraction intensities.

# the conventional unit cell is given here
indexed_spots = measurement.index_diffraction_spots(cell=sio2.cell)

We show the last diffraction spots in the thickness series with an overlay showing the miller indices.

indexed_spots[-1].block_direct().show(scale=0.5, figsize=(5, 5));
../../../_images/f62fd9b9f32d3c1b93e54e99fb7b0db927e748707625b4291f698b476b86fbd6.png

The diffraction spots can be converted to a pandas dataframe.

import pandas as pd

pd.set_option("display.max_columns", 8)

df = indexed_spots.to_dataframe()

df
-5 -6 0 -5 -5 0 -5 -4 0 -5 -3 0 ... 5 3 0 5 4 0 5 5 0 5 6 0
0.000000 0.000000 0.000000 0.000000 0.000000 ... 0.000000 0.000000 0.000000 0.000000
99.972725 0.000532 0.001930 0.000604 0.001943 ... 0.001902 0.000548 0.001845 0.000489
199.945451 0.001056 0.012364 0.002705 0.001561 ... 0.001537 0.002574 0.012001 0.000941
299.918176 0.001704 0.016678 0.003978 0.000847 ... 0.000844 0.003900 0.016278 0.001630
399.890901 0.001201 0.010127 0.002770 0.000779 ... 0.000736 0.002690 0.009888 0.001215

5 rows × 120 columns

visualization = indexed_spots[-1].block_direct().show(scale=0.5)
visualization.axis_off(spines=False)
visualization.fig.savefig(
    "../thumbnails/saed_quickstart.png", bbox_inches="tight", pad_inches=0
)
../../../_images/4766001ad348c0a28ae704ff7407fe5d8db6c8d43d39a8308736e07787c07d5a.png