Analysing atom lattices

After finding and refining the atom lattices as shown in Finding the atom lattice, the atomic structure can be analysed through

  1. Ellipticity of the atomic columns

  2. Monolayer separation

  3. Angle between monolayers

  4. Integration of atomic columns

  5. Making line profiles

  6. Plotting a pair distribution function

In this tutorial we will use a dummy image containing two sublattices. Different structural distortions have been introduced in the image, and this tutorial will study these distortions. The plots that are made are simple plots that clearly shows the structural changes in the sublattices. It is also possible to make more fancy figures for publications (such as in the paper presenting Atomap), and some tips for doing this are shown in the plotting section of the guide.

Fantasite, the dummy structure

The procedure for finding and refining the sublattices in fantasite is similar as in Images with more than one sublattice.

>>> import atomap.api as am
>>> from atomap.tools import remove_atoms_from_image_using_2d_gaussian

>>> s = am.dummy_data.get_fantasite()
>>> A_positions = am.get_atom_positions(s, separation=12, pca=True)
>>> sublattice_A = am.Sublattice(A_positions, image=s.data, color='r', name='A')
>>> sublattice_A.find_nearest_neighbors()
>>> sublattice_A.refine_atom_positions_using_center_of_mass()
>>> sublattice_A.refine_atom_positions_using_2d_gaussian()
>>> sublattice_A.construct_zone_axes()

>>> direction_001 = sublattice_A.zones_axis_average_distances[1]
>>> B_positions = sublattice_A.find_missing_atoms_from_zone_vector(direction_001)
>>> image_without_A = remove_atoms_from_image_using_2d_gaussian(sublattice_A.image, sublattice_A)

>>> sublattice_B = am.Sublattice(B_positions, image_without_A, color='blue', name='B')
>>> sublattice_B.construct_zone_axes()
>>> sublattice_B.refine_atom_positions_using_center_of_mass()
>>> sublattice_B.refine_atom_positions_using_2d_gaussian()
>>> atom_lattice = am.Atom_Lattice(image=s.data, name='fantasite', sublattice_list=[sublattice_A, sublattice_B])
>>> atom_lattice.save("fantasite.hdf5", overwrite=True)
>>> s.plot()
>>> atom_lattice.plot()
_images/fantasite.png _images/atom_lattice.png

Fantasite is shown in the left image, and it is possible to see some variations in ellipticity and atom positions with the naked eye. The right figure shows the atom positions after the refinement. As refinement of the sublattice can be time consuming, it is a good idea to save the final atom lattice.

The Atom_Lattice object

The atom lattice can be loaded:

>>> import atomap.api as am
>>> atom_lattice = am.load_atom_lattice_from_hdf5("fantasite.hdf5")
>>> atom_lattice
<Atom_Lattice, fantasite (sublattice(s): 2)>
>>> atom_lattice.sublattice_list 
[<Sublattice, A (atoms:497,planes:6)>, <Sublattice, B (atoms:465,planes:6)>] 
>>> image = atom_lattice.image

atomap.atom_lattice.Atom_Lattice is an object containing the sublattices, and other types of information. The fantasite atom lattice contains two sublattices (red and blue dots in the image above). Atom positions, sigma, ellipticity and rotation for the atomic columns in a sublattice can be accessed via attributes in the sublattice object.

>>> sublattice_A = atom_lattice.sublattice_list[0]
>>> x = sublattice_A.x_position
>>> y = sublattice_A.y_position
>>> sigma_x = sublattice_A.sigma_x
>>> sigmal_y = sublattice_A.sigma_y
>>> ellipticity = sublattice_A.ellipticity
>>> rotation = sublattice_A.rotation_ellipticity

Similarly, properties of a single atomic column atomap.atom_position.Atom_Position can be accessed in atomap.sublattice.Sublattice.atom_list. The atomap.atom_position.Atom_Position object contain information related to a specific atomic column.

>>> atom_position_list = sublattice_A.atom_list
>>> atom_position = atom_position_list[0]
>>> x = atom_position.pixel_x
>>> y = atom_position.pixel_y
>>> sigma_x = atom_position.sigma_x
>>> sigma_y = atom_position.sigma_y

The atomap.atom_plane.Atom_Plane objects contain the atomic columns belonging to the same specific plane. Atom plane objects are defined by the direction vector parallel to the atoms in the plane, for example (58.81, -41.99). These can be accessed by:

>>> atom_plane_list = sublattice_A.atom_plane_list
>>> atom_plane = atom_plane_list[0]
>>> atoms_in_plane_list = atom_plane.atom_list

Atom planes for a specific direction can also be accessed:

>>> zone_axis0 = sublattice_A.zones_axis_average_distances[0]
>>> atom_plane_list0 = sublattice_A.atom_planes_by_zone_vector[zone_axis0]

Ellipticity

Elliptical atomic columns may occur when atoms parallel to the electron beam have shifted position in the plane perpendicular to the beam. In the image, circular atomic columns have an ellipticity (\epsilon) of 1, as sigma_x = sigma_y (\sigma_x = \sigma_y). Ellipticity is defined as

\epsilon =
    \begin{cases}
            \frac{\sigma_x}{\sigma_y},& \text{if } \sigma_x > \sigma_y\\
                    \frac{\sigma_y}{\sigma_x},& \text{if } \sigma_y > \sigma_x\\
                        \end{cases}

Ellipticity maps

The ellipticity map shows the magnitude of the ellipticity. Values are interpolated, giving a continuous map. The sublattice B was generated without any ellipticity, and the image to the right is fairly flat, as expected. In sublattice A, a region with elliptical atomic columns is clearly visible. The ellipticity gradient is perfectly in line with how the dummy image of fantasite has been generated. Maps gives nice visualization of gradual change.

>>> sublattice_A = atom_lattice.sublattice_list[0]
>>> sublattice_B = atom_lattice.sublattice_list[1]
>>> sublattice_A.plot_ellipticity_map(cmap='viridis', vmin=0.95, vmax=1.3)
>>> sublattice_B.plot_ellipticity_map(cmap='viridis', vmin=0.95, vmax=1.3)
_images/ellipticity_map_A.png _images/ellipticity_map_B.png

The atomap.sublattice.Sublattice.plot_ellipticity_map() function calls atomap.sublattice.Sublattice.get_ellipticity_map(), which calls atomap.sublattice.Sublattice.get_property_map().

Vector plots

While the ellipticity map nicely visualizes the magnitude (and gradual change) of the ellipticity, it does not show the direction of the ellipticity. In vector (quiver) plots (plot_ellipticity_vectors()) both the rotation and magnitude are visualized, through the length and angle of the arrows. There is one arrow for each atom position. This type of visualization also reveals that the ellipticity of the atomic columns have an alternating directionality.

>>> sublattice_A.plot_ellipticity_vectors()
_images/ellipticity_vectors.png

In this function, a value of 1 is subtracted from the magnitude of the ellipticity. This makes it easier to study changes in ellipticity, as the 0-point of the plot is set to the perfect circle. atomap.plotting.plot_vector_field() is called, and this function uses Matplotlib’s quiver plot function.

Distance between monolayers

As Atomap knows the positions of all atoms, it can also tell you if you have any variations in the unit cell sizes or other types of structural distortions. For example, Atomap has been used to study oxygen octahedron tilt patterns in perovskite thin films.

In this example, Atomap finds the distance between monolayers. “Distance between monolayers” is defined in the figure below (a). The distance is measured between an atomic column and the nearest monolayer, as shown in this figure. The values for monolayer separation is attributed to the midpoint of the atom and monolayer.

_images/definition.jpg

This is found by using the get_monolayer_distance_map(), which returns a stack of HyperSpy images. In this signal the navigation axis is a zone vector (the values in the zone axis tuple are in principle the vector parallel to the monolayer) and the signal axes shows monolayer separation at each atom position. get_monolayer_distance_map() can also take in a subset of zone vectors, but the default is to find the monolayer separation for all the zone axes.

>>> s_monolayer = sublattice_B.get_monolayer_distance_map()
>>> s_monolayer.plot(cmap='viridis')
_images/Sublattice_B_monolayer_distance_a.png _images/Sublattice_B_monolayer_distance_b.png

The left image shows the monolayer separation for one zone axis, namely the separation between the monolayers drawn up by red lines in the right figure. Clearly, the position of the B atomic columns are changed in the middle of the image, where every second monolayer is closer and further apart from the atom.

Atom distance difference

An alternative to plotting the monolayer separation is the difference in distance from an atom to the two neighbouring atoms in the monolayer with get_atom_distance_difference_map(). In some cases, this is a more robust way of finding superstructures as this method only shows the variations in distance, not the total distance. The “distance difference” is defined in the below figure. The value of the distance difference is attributed to the position of the middle atom.

_images/definition2.jpg

The B-sublattice of fantasite has two different domains of structural distortion. As seen above, the monolayers in (0.0, -30.0) direction are alternatively closer and further apart. The corresponding distance difference will be the distance difference for atoms in the orthogonal planes (index 0, value (15.0, 0.0), shown below).

>>> zone = sublattice_B.zones_axis_average_distances[0]
>>> s_dd = sublattice_B.get_atom_distance_difference_map([zone])
>>> s_dd.plot(cmap='viridis')
>>> sublattice_B.plot_planes(image=atom_lattice.image)
_images/Angle_map_z1.png _images/sublatticeB_dd_map_0.png

The left image shows the monolayers, with the resulting distance difference maps shown in the image to the right.

The second structural domain is found in the vertical planes.

>>> zone = sublattice_B.zones_axis_average_distances[1]
>>> s_dd = sublattice_B.get_atom_distance_difference_map([zone], add_zero_value_sublattice=sublattice_A)
>>> s_dd.plot(cmap='viridis')
_images/Angle_map_z2.png _images/sublatticeB_dd_map_1.png

Angle between atoms

To visualize the angle between the atoms, get_atom_angles_from_zone_vector() is used. This function returns three lists: x- and y-coordinates of the atoms, and a list of the angle between two zone axes at each atom.

>>> z1 = sublattice_B.zones_axis_average_distances[0]
>>> z2 = sublattice_B.zones_axis_average_distances[1]
>>> x, y, a = sublattice_B.get_atom_angles_from_zone_vector(z1, z2, degrees=True)
>>> s_angle = sublattice_B.get_property_map(x, y, a)
>>> s_angle.plot()
_images/Angle_map.png _images/Angle_map_zoom.png

The atomic columns start to “zigzag” in the rightmost part of the image. This is also clear with the naked eye (atomic columns marked with blue dots). get_property_map() is a very general function, and can plot a map of any property.

Integration of atomic columns

When analysing the intensity of different atomic columns it is important to be able to accurately integrated over all columns in the field of view in an automated way. Two methods of image segmentation have been implemented into Atomap, these are Voronoi cell integration and Watershed. The Voronoi method can be applied to 3D data-sets e.g. EDX and EELS.

The integrate() function returns a list containing:

  1. Integrated intensity list - the same length as atom coordinates list with resulting integrated intensities.

  2. Intensity record - Image the same size as input image, each pixel in a particular segment or region has the value of the integration of the region.

  3. Point record - An image where defining the locations of each integration region.

>>> i_points, i_record, p_record = atom_lattice.integrate_column_intensity()
>>> i_record.plot()
_images/intensity_record_voronoi.png

To use the Watershed method:

>>> i_points, i_record, p_record = atom_lattice.integrate_column_intensity(
...         method='Watershed')
>>> i_record.plot()
_images/intensity_record_watershed.png

The Voronoi cell integration has a max_radius optional input which helps to prevent cells from becoming too large e.g. at the surface of a nanoparticle.

Line profiles

Often it can be a good idea to integrate parts of the image, for example to improve the signal-to-noise ratio or reduce the information to fewer dimensions. This example will introduce how line profiles of properties can be made.

>>> zone = sublattice_A.zones_axis_average_distances[1]
>>> plane = sublattice_A.atom_planes_by_zone_vector[zone][23]
>>> s_elli_line = sublattice_A.get_ellipticity_line_profile(plane)
_images/line_ellip.png _images/line_ellip_plane.png

In the above example, plane is atom plane number 23 in the right figure. This the “last plane” from the left where the atomic columns in sublattice_A are circular. In plane 22, the columns become elliptical. This plane is therefore sent into the function get_ellipticity_line_profile(), as the definition of where the interface (between two structural domains) is. The values for ellipticity are integrated and averaged in the same direction as the plane (downwards), and plotted in the left figure. Average ellipticity is on the y-axis, while the x-axis is the distance from the interface (plane 23). On the x-axis, negative values are to the left of the interface and positive values to the right of the interface.

get_ellipticity_line_profile() calls _get_property_line_profile(), which takes in 3 lists: x and y coordinates for the atoms, and a list of value for a property (in this case ellipticity). It then sorts the atoms after distance from interface. The atoms with the same distance from the interface belong to the same plane, parallel to the interface, and the value for a property (such as ellipticity) for these atoms are averaged.

The raw data can be accessed in the signal’s metadata:

>>> position = s_elli_line.metadata.line_profile_data.x_list
>>> ellipticity = s_elli_line.metadata.line_profile_data.y_list
>>> standard_deviation = s_elli_line.metadata.line_profile_data.std_list

And plotted using matplotlib, with the standard deviation shown in red.

>>> import matplotlib.pyplot as plt
>>> fig, ax = plt.subplots()
>>> _ = ax.errorbar(position, ellipticity, yerr=standard_deviation, ecolor='red')
>>> fig.show()
_images/line_profile_errorbar.png

There are also functions to make such line profiles for monolayer separation (get_monolayer_distance_line_profile()), and “distance difference” (get_atom_distance_difference_line_profile()). These two properties are “directional”, which means that the zone axis for the distance measurement must be given in addition to the “interface” plane.

>>> zone = sublattice_B.zones_axis_average_distances[1]
>>> plane = sublattice_B.atom_planes_by_zone_vector[zone][0]
>>> s_monolayer_line = sublattice_B.get_monolayer_distance_line_profile(zone, plane)
_images/line_monolayer.png

The property of interest here was the separation of monolayers, and the separation between two specific monolayers is plotted at the point between these two monolayers.

>>> zone = sublattice_B.zones_axis_average_distances[1]
>>> plane = sublattice_B.atom_planes_by_zone_vector[zone][-1]
>>> s_dd_line = sublattice_B.get_atom_distance_difference_line_profile(zone, plane)
_images/line_dd.png

Naturally, in these functions any atomic plane can be defined as the interface.

Finding polarization

In many ferroelectric materials, the spontaneous electric polarization can be determined by looking at the shift of some atomic columns in relation to the others. One example of this is in the ferroelectric perovskite oxides, where the B-cation is shifted from its cubic centrosymmetric position. The polarization can be determined by finding both the direction and magnitude of this shift.

Firstly, we get an appropriate artificial dataset, resembling a ferroelectric thin film grown on top of a non-ferroelectric substrate.

>>> atom_lattice = am.dummy_data.get_polarization_film_atom_lattice()
>>> atom_lattice.plot()
_images/polarization_atom_lattice.png

The blue, B-cation, atom columns in the top part of the image are shifted towards the left. By finding the centre position of four neighboring red A-cation forming a square, this shift can be quantified.

Finding these neighbors relies on moving along two zone axis directions in the A-cation sublattice.

>>> sublatticeA = atom_lattice.sublattice_list[0]
>>> sublatticeA.construct_zone_axes()

Next, find the two perpendicular zone axes spanning this square. For the perovskite oxide 100 projection, this is most likely the two first ones.

>>> sublatticeA.plot_planes()
_images/polarization_atom_plane0.png _images/polarization_atom_plane1.png

The zone axes are then used with the B-cation sublattice in the get_polarization_from_second_sublattice() method:

>>> za0 = sublatticeA.zones_axis_average_distances[0]
>>> za1 = sublatticeA.zones_axis_average_distances[1]
>>> sublatticeB = atom_lattice.sublattice_list[1]
>>> s_polarization = sublatticeA.get_polarization_from_second_sublattice(za0, za1, sublatticeB)

This can be visualized directly by using the plot method, and the data itself can be accessed in the signal’s metadata.

>>> s_polarization.plot()
>>> vector_list = s_polarization.metadata.vector_list
_images/polarization_signal_marker.png

Plotting a pair distribution function

A method of analysis common in X-ray diffraction is the plotting of a pair distribution function (PDF). This function calculates the distance between every atom with every other atom in the sampled volume. This function can also be performed in two dimensions from an image of atoms. A PDF of a crystalline material is similar to a radially averaged diffraction pattern. They are particularly useful for looking for ordering in amorphous materials.

A PDF can be calculated from a Sublattice object using pair_distribution_function():

>>> sublattice = am.dummy_data.get_simple_cubic_sublattice()
>>> s_pdf = sublattice.pair_distribution_function()
>>> s_pdf.plot()
_images/pair_distribution_function.png