Show code cell source
%config InlineBackend.rc = {"figure.dpi": 72, "figure.figsize": (6.0, 4.0)}
%matplotlib inline
import ase
import matplotlib.pyplot as plt
import numpy as np
from ase.build import bulk
import abtem
CBED quickstart#
This notebook demonstrates a basic CBED simulation of silicon in the \((111)\) zone axis.
Configuration#
We start by (optionally) setting our configuration. See documentation for details.
abtem.config.set({"device": "cpu", "fft": "fftw"})
<abtem.core.config.set at 0x22d18f4ff70>
Atomic model#
We create the atomic model. See our walkthough or our tutorial on atomic models.
We create an atomic model of silicon using the bulk
function from ASE.
silicon = bulk("Si", crystalstructure="diamond")
abtem.show_atoms(silicon, plane="xy");
![../../../_images/91a8c48d713a6f62d00c9b8f85eab035861192b2682db1ca56f4388662796fb8.png](../../../_images/91a8c48d713a6f62d00c9b8f85eab035861192b2682db1ca56f4388662796fb8.png)
We can choose the \((111)\) zone axis to be a propagation direction using the surface
function.
silicon_111 = ase.build.surface(
silicon, (1, 1, 1), layers=3, periodic=True
) # create surface structure in the (111) direction
silicon_111_orthogonal = abtem.orthogonalize_cell(silicon_111) # make cell orthogonal
abtem.show_atoms(silicon_111_orthogonal);
![../../../_images/fff944d3eba2f7b11eac4bb12c59e58da44765fc7b2f44a3c3d8a8209a4c370f.png](../../../_images/fff944d3eba2f7b11eac4bb12c59e58da44765fc7b2f44a3c3d8a8209a4c370f.png)
Finally we repeat the structure in \(x\) and \(y\), this will improve the reciprocal space resolution. We also repeat the structure in \(z\), thus simulating a thicker sample.
atoms = silicon_111_orthogonal * (13, 8, 60)
fig, (ax1, ax2) = plt.subplots(1, 2, figsize=(12, 6))
abtem.show_atoms(atoms, ax=ax1, title="Beam view")
abtem.show_atoms(atoms, ax=ax2, plane="xz", title="Side view", linewidth=0.0);
![../../../_images/ae8049eebde6e91520832a05969a378eae584ce3ecbe8cfa6b91eee3f2318aa8.png](../../../_images/ae8049eebde6e91520832a05969a378eae584ce3ecbe8cfa6b91eee3f2318aa8.png)
Potential#
We create an ensemble of potentials using the frozen phonon. See our walkthrough on frozen phonons.
frozen_phonons = abtem.FrozenPhonons(atoms, 8, {"Si": 0.078})
We create a potential from the frozen phonons model, see walkthrough on potentials.
potential = abtem.Potential(
frozen_phonons,
sampling=0.2,
projection="infinite",
slice_thickness=2,
exit_planes=80,
)
Wave function#
We create a probe wave function at an energy of 100 keV with a convergence semiangle of \(9.4 \ \mathrm{mrad}\). See our walkthrough on wave functions.
Partial temporal coherence is neglected here, see our tutorial on partial coherence if you want to include this in your simulation..
wave = abtem.Probe(energy=100e3, semiangle_cutoff=9.4)
wave.grid.match(potential)
wave.profiles().show();
[########################################] | 100% Completed | 224.81 ms
![../../../_images/7b44cd4af068be51a6f4c186a46ba8954889c5856e3b3738e4567000ccb37e9b.png](../../../_images/7b44cd4af068be51a6f4c186a46ba8954889c5856e3b3738e4567000ccb37e9b.png)
Multislice#
We run the multislice algorithm and calculate the diffraction patterns, see our walkthrough on multislice.
measurement = wave.multislice(potential).diffraction_patterns(max_angle=30)
We take the mean across the frozen phonons axis, and compute the result.
measurement = measurement.mean(0)
measurement.compute()
[########################################] | 100% Completed | 11.78 ss
<abtem.measurements.DiffractionPatterns object at 0x0000022D1D471090>
Visualize results#
We show the thickness series as an exploded plot.
visualization = measurement.show(
explode=True,
figsize=(16, 5),
)
![../../../_images/1ec5a29cea2831c9b43ffd0fb49ec39abe5279a5d09aaa1d965af5f0e65842c6.png](../../../_images/1ec5a29cea2831c9b43ffd0fb49ec39abe5279a5d09aaa1d965af5f0e65842c6.png)