API documentation
Classes
Atom lattice
Dumbbell lattice
Sublattice
Atom plane
- class atomap.atom_plane.Atom_Plane(atom_list, zone_vector, atom_lattice)
- Parameters:
- property amplitude_gaussian
- property amplitude_max_intensity
- property ellipticity
- get_angle_to_horizontal_axis()
Get angle between atoms in the atom plane and horizontal axis.
- get_atom_distance_list()
- get_atom_index(check_atom)
- get_closest_distance_and_angle_to_point(points_x, points_y, use_precalculated_line=False, plot_debug=False)
Return the smallest distances from each point in a list to the atom plane.
- Parameters:
points_x (list of numbers) – The x and y coordinates for the atoms.
points_y (list of numbers) – The x and y coordinates for the atoms.
use_precalculated_line (bool or list/array) – The coefficients [a, b] for the linear line y = ax + b to which the shortest distance should be found. By default False, in which the coefficients are found by fitting a straight line to self (atom plane).
plot_debug (bool, default False) – If true, a debug plot is saved. The plot is a 3D plot of x and y positions and distance to the plane.
- Returns:
list of numbers
- Return type:
list of the shortest distances.
- get_closest_position_to_point(point_position, extend_line=False)
- get_connecting_atom_planes(atom_plane, zone_vector)
Get the outer atom planes which connect self atom plane with another atom plane, through a specific atom plane direction.
The self atom plane, atom plane, and the two returned atom planes will span make a parallelogram.
- get_intersecting_atom_from_atom_plane(atom_plane)
- get_net_distance_change_between_atoms()
Output [(x,y,z)]
- get_slice_between_two_atoms(atom1, atom2)
- get_x_position_list()
- get_y_position_list()
- position_distance_to_neighbor()
Get distance between all atoms and its next neighbor in the atom plane, and the position between these two atoms.
- Return type:
Numpy array [x, y, distance]
Example
>>> from numpy.random import random >>> from atomap.sublattice import Sublattice >>> pos = [[x, y] for x in range(9) for y in range(9)] >>> sublattice = Sublattice(pos, random((9, 9))) >>> sublattice.construct_zone_axes() >>> atom_plane = sublattice.atom_plane_list[10] >>> pos_distance = atom_plane.position_distance_to_neighbor() >>> x_pos = pos_distance[:, 0] >>> y_pos = pos_distance[:, 1] >>> distance = pos_distance[:, 2]
- property rotation
- property rotation_ellipticity
- property sigma_average
- property sigma_x
- property sigma_y
- sort_atoms_by_distance_to_point(point=(0, 0))
- property x_position
- property y_position
Atom position
- class atomap.atom_position.Atom_Position(x, y, sigma_x=1.0, sigma_y=1.0, rotation=0.01, amplitude=1.0)
The Atom_Position class contain information about a single atom column.
- Parameters:
x (float)
y (float)
sigma_x (float, optional)
sigma_y (float, optional)
rotation (float, optional) – In radians. The rotation of the axes of the 2D-Gaussian relative to the image axes. In other words: the rotation of the sigma_x relative to the horizontal axis (x-axis). This is different from the rotation_ellipticity, which is the rotation of the largest sigma in relation to the horizontal axis. For the rotation of the ellipticity, see rotation_ellipticity.
amplitude (float, optional) – Amplitude of Gaussian. Stored as amplitude_gaussian attribute.
- rotation
The rotation of sigma_x axis, relative to the x-axis in radians. This value will always be between 0 and pi, as the elliptical 2D-Gaussian used here is always symmetric in the rotation direction, and the perpendicular rotation direction.
- Type:
- rotation_ellipticity
The rotation of the longest sigma, relative to the x-axis in radians.
- Type:
- refine_position
If True (default), the atom position will be fitted to the image data when calling the Sublattice.refine_atom_positions_using_center_of_mass and Sublattice.refine_atom_positions_using_2d_gaussian methods. Note, the atom will still be fitted when directly calling the Atom_Position.refine_position_using_center_of_mass and Atom_Position.refine_position_using_2d_gaussian methods in the atom position class itself. Setting it to False can be useful when dealing with vacancies, or other features where the automatic fitting doesn’t work.
- Type:
Examples
>>> from atomap.atom_position import Atom_Position >>> atom_position = Atom_Position(10, 5)
More parameters
>>> atom_pos = Atom_Position(10, 5, sigma_x=2, sigma_y=4, rotation=2)
- as_gaussian()
- calculate_max_intensity(image_data, percent_to_nn=0.4)
Find the maximum intensity of the atom. See get_atom_column_amplitude_max_intensity() for further uses.
The maximum intensity is found within the distance to the nearest neighbor times percent_to_nn.
- Parameters:
image_data (NumPy 2D array)
percent_to_nn (float, default 0.4) – Determines the boundary of the area surrounding each atomic column, as fraction of the distance to the nearest neighbour.
- Return type:
Maximum pixel intensity for an atom position.
Example
>>> import atomap.api as am >>> sublattice = am.dummy_data.get_simple_cubic_sublattice() >>> sublattice.find_nearest_neighbors() >>> atom0 = sublattice.atom_list[0] >>> atom0_max_int = atom0.calculate_max_intensity(sublattice.image)
- calculate_min_intensity(image_data, percent_to_nn=0.4)
Find the minimum intensity of the atom. See get_atom_column_amplitude_min_intensity() for further uses.
The min intensity is found within the distance to the nearest neighbor times percent_to_nn.
- Parameters:
image_data (NumPy 2D array)
percent_to_nn (float, default 0.4) – Determines the boundary of the area surrounding each atomic column, as fraction of the distance to the nearest neighbour.
- Return type:
Minimum pixel intensity for an atom position.
Example
>>> import atomap.api as am >>> sublattice = am.dummy_data.get_simple_cubic_sublattice() >>> sublattice.find_nearest_neighbors() >>> atom0 = sublattice.atom_list[0] >>> atom0_min_int = atom0.calculate_min_intensity(sublattice.image)
- can_atom_plane_be_reached_through_zone_vector(atom_plane, zone_vector)
- property ellipticity
Largest sigma divided by the shortest
- estimate_local_scanning_distortion(image_data, radius=6, edge_skip=2)
Get the amount of local scanning distortion from an atomic column.
This is done by assuming the atomic column has a symmetrical shape, like Gaussian or Lorentzian. The distortion is calculated by getting the image_data around the atom, given by the radius parameter. This gives a square cropped version of the image_data, where the region outside the radius is masked. For each line in the horizontal direction, the center of mass is found. This gives a list of horizontal positions as a function of the vertical lines. To remove the effects like astigmatism and mistilt a linear fit is fitted to this list of horizontal positions. This fit is then subtracted from the horizontal position. The distortion for the vertical lines is then calculated by getting the standard deviation of this list of values. Getting the horizontal distortions is calculated in a similar fashion, but with a list of vertical positions as a function of the horizontal lines.
- Parameters:
image_data (2D NumPy array)
radius (int) – Radius of the masked and cropped image. Default 6.
edge_skip (int) – When the cropped image is masked with a circle, the edges will consist of very few unmasked pixels. The center of mass of these lines will be unreliable, so they’re by default skipped. edge_skip = 2 means the two lines closest to the edge of the cropped image will be skipped. Default 2.
- Returns:
scanning_distortion – Both horizontal and vertical directions. For standard raster scans, horizontal will be the fast scan direction, and y the slow scan direction. Thus the returned values will be: (fast scan, slow scan), where typically the slow scan direction has the largest amount of distortion.
- Return type:
Examples
>>> sl = am.dummy_data.get_scanning_distortion_sublattice() >>> atom = sl.atom_list[50] >>> distortion = atom.estimate_local_scanning_distortion(sl.image)
- find_atom_intensity_inside_mask(image_data, radius)
Find the average intensity inside a circle.
The circle is defined by the atom position, and the given radius (in pixels). The outside this area is covered by a mask. The average intensity is saved to self.intensity_mask.
- get_angle_between_atoms(atom0, atom1=None)
Return the angle between atoms in radians.
Can either find the angle between self and two other atoms, or between another atom and the horizontal axis.
- Parameters:
atom0 (Atom Position object) – The first atom.
atom1 (Atom Position object, optional) – If atom1 is not specified, the angle between itself, atom0 and the horizontal axis will be returned.
- Returns:
Angle – Angle in radians
- Return type:
Examples
>>> from atomap.atom_position import Atom_Position >>> atom0 = Atom_Position(0, 0) >>> atom1 = Atom_Position(1, 1) >>> atom2 = Atom_Position(-1, 1) >>> angle0 = atom0.get_angle_between_atoms(atom1, atom2) >>> angle1 = atom0.get_angle_between_atoms(atom1)
- get_angle_between_zone_vectors(zone_vector0, zone_vector1)
Return the angle between itself and the next atoms in the atom planes belonging to zone_vector0 and zone_vector1
- get_atomic_plane_from_zone_vector(zone_vector)
- get_closest_neighbor()
Find the closest neighbor to an atom in the same sub lattice.
- Return type:
Atomap atom_position object
- get_ellipticity_rotation_vector()
- get_ellipticity_vector()
- get_index_in_atom_plane(atom_plane)
- get_neighbor_atoms_in_atomic_plane_from_zone_vector(zone_vector)
- get_next_atom_in_atom_plane(atom_plane)
- get_next_atom_in_zone_vector(zone_vector)
Get the next atom in the atom plane belonging to zone vector.
- get_pixel_difference(atom)
Vector between self and given atom
- get_pixel_distance_from_another_atom(atom)
- get_pixel_position()
- get_position_convergence(distance_to_first_position=False)
- get_previous_atom_in_atom_plane(atom_plane)
- get_previous_atom_in_zone_vector(zone_vector)
- get_rotation_vector()
- is_in_atomic_plane(zone_direction)
- pixel_distance_from_point(point=(0, 0))
- refine_position_using_2d_gaussian(image_data, rotation_enabled=True, percent_to_nn=0.4, mask_radius=None, centre_free=True)
Use 2D Gaussian to refine the parameters of the atom position.
- Parameters:
image_data (Numpy 2D array)
rotation_enabled (bool, optional) – If True, the Gaussian will be able to rotate. Note, this can increase the chance of fitting failure. Default True.
percent_to_nn (float, optional) – The percent of the distance to the nearest neighbor atom in the same sub lattice. The distance times this percentage defines the mask around the atom where the Gaussian will be fitted. A smaller value can reduce the effect from neighboring atoms, but might also decrease the accuracy of the fitting due to less data to fit to. Default 0.4 (40%).
mask_radius (float, optional) – Radius of the mask around each atom. If this is not set, the radius will be the distance to the nearest atom in the same sublattice times the percent_to_nn value. Note: if mask_radius is not specified, the Atom_Position objects must have a populated nearest_neighbor_list.
centre_free (bool, default True) – If True, the centre parameter will be free, meaning that the Gaussian can move.
- refine_position_using_center_of_mass(image_data, percent_to_nn=0.4, mask_radius=None)
Refine the position of the atom position using center of mass
The position is stored in atom_position.pixel_x and atom_position.pixel_y. The old positions are saved in atom_position.old_pixel_x_list and atom_position.old_pixel_x_list.
- Parameters:
image_data (Numpy 2D array)
percent_to_nn (float, optional) – The percent of the distance to the nearest neighbor atom in the same sub lattice. The distance times this percentage defines the mask around the atom where the Gaussian will be fitted. A smaller value can reduce the effect from neighboring atoms, but might also decrease the accuracy of the fitting due to less data to fit to. Default 0.4 (40%).
mask_radius (float, optional) – Radius of the mask around each atom. If this is not set, the radius will be the distance to the nearest atom in the same sublattice times the percent_to_nn value. Note: if mask_radius is not specified, the Atom_Position objects must have a populated nearest_neighbor_list.
Examples
>>> from atomap.atom_position import Atom_Position >>> atom = Atom_Position(15, 10) >>> image_data = np.random.randint(100, size=(20, 20)) >>> atom.refine_position_using_center_of_mass( ... image_data, mask_radius=5)
- property rotation
The rotation of the atom relative to the horizontal axis.
Given in radians. For the rotation of the ellipticity, see rotation_ellipticity.
- property rotation_ellipticity
Rotation between the “x-axis” and the major axis of the ellipse.
Rotation between the horizontal axis, and the longest part of the atom position, given by the longest sigma. Basically giving the direction of the ellipticity.
- set_element_info(element, z)
Set which atoms are present along the atomic column at this atom position.
- Parameters:
Examples
>>> from atomap.atom_position import Atom_Position >>> atom = Atom_Position(x=15, y=10, sigma_x=5, sigma_y=3) >>> atom.set_element_info("C", [0, 2.5]) >>> atom2 = Atom_Position(x=10, y=15, sigma_x=5, sigma_y=3) >>> atom2.set_element_info(["C", "O"], [0, 2.5])
- property sigma_average
- property sigma_x
- property sigma_y
Process parameters
- class atomap.process_parameters.ModelParametersBase
- add_sublattice_config(sublattice_config_object)
- get_sublattice_from_order(order_number)
- property number_of_sublattices
- class atomap.process_parameters.GenericStructure
A Generic structure with one sublattice, the GenericSublattice.
- Parameters:
peak_separation (None) – No peak separation is set, pixel_separation must be given separately in order to find initial atom positions.
name ('A structure')
sublattice_list (list of Sublattices) – Contains one GenericSublattice
- class atomap.process_parameters.PerovskiteOxide110
The Perovskite Oxide structure in the 110 projection.
- Parameters:
name ('Perovskite 110')
peak_separation (0.127 , distance in nm) – Approximately half the distance between the most intense atoms in the structure, used to get initial atom position by peak finding.
sublattice_list (list of Sublattices) – Contains 3 sublattices, for A, B and O: PerovskiteOxide110SublatticeACation, PerovskiteOxide110SublatticeBCation, PerovskiteOxide110SublatticeOxygen
- class atomap.process_parameters.SrTiO3_110
- class atomap.process_parameters.SublatticeParameterBase
- class atomap.process_parameters.GenericSublattice
Process parameters for a generic sublattice
- color
- Type:
‘red’ , color of the markers indicating atom positions
- image_type
The image will not be inverted.
- Type:
0
- name
- Type:
‘S0’ , name of the sublattice
- sublattice_order
The sublattice will get the order 0. Higher orders can be used for less intense sublattices in images with multiple sublattices.
- Type:
0
- refinement_config
Dict with configuration settings for the refinement of atom positions. 1st refinement : center-of-mass on image modified with background removal, PCA noise filtering and normalization. 2nd refinement : Center of mass on the original image. 3rd refinement : Fitting 2D Gaussians to the original image.
- Type:
- neighbor_distance
Mask radius for fitting, set to 0.35 of the distance to nearest neighbor.
- Type:
0.35
- atom_subtract_config
Configuration for how to subtract higher order sublattices from the image. (Not really used in this structure, but included as an example.)
- Type:
Examples
>>> import atomap.process_parameters as pp >>> generic_sublattice = pp.GenericSublattice()
- class atomap.process_parameters.PerovskiteOxide110SublatticeACation
Process parameters for the most intense atoms (typically the A-cations) in a Perovskite Oxide structure projected along the (110) direction.
- color
- Type:
‘blue’ , color of the markers indicating atom positions
- image_type
The image will not be inverted for finding atom positions.
- Type:
0
- name
- Type:
‘A-cation’ , name of the sublattice
- sublattice_order
The most intense sublattice gets order 0.
- Type:
0
- refinement_config
Dict with configuration settings for the refinement of atom positions. Two refinements by fitting 2D Gaussians to the original image.
- Type:
- neighbor_distance
Mask radius for fitting set to 0.35 of nearest neighbor.
- Type:
0.35
- class atomap.process_parameters.PerovskiteOxide110SublatticeBCation
Process parameters for the second most intense atoms (typically the B-cations) in a Perovskite Oxide structure projected along the (110) direction.
- color
- Type:
‘green’ , color of the markers indicating atom positions
- image_type
The image will not be inverted.
- Type:
0
- name
- Type:
‘B-cation’ , name of the sublattice
- sublattice_order
The sublattice will get the order 1. Higher orders can be used for less intense sublattices in images with multiple sublattices. Lower order is for more intense sublattices.
- Type:
1
- refinement_config
Dict with configuration settings for the refinement of atom positions. 1st refinement : center-of-mass on the original image. 2nd refinement : Fitting 2D Gaussians to the original image.
- Type:
- neighbor_distance
Mask radius for fitting set to 0.25 of nearest neighbor.
- Type:
0.25
- class atomap.process_parameters.PerovskiteOxide110SublatticeOxygen
Process parameters for the least intense atoms (typically the Oxygen) in a Perovskite Oxide structure projected along the (110) direction.
- color
- Type:
‘red’ , color of the markers indicating atom positions
- image_type
The image will be inverted. Oxygen can be visible in ABF images. The image is inverted such that the atoms are bright dots on a dark background.
- Type:
1
- name
- Type:
‘Oxygen’ , name of the sublattice
- sublattice_order
The sublattice will get the order 2, being the third most intense sublattice.
- Type:
2
- refinement_config
Dict with configuration settings for the refinement of atom positions. 1st refinement : center-of-mass on the original image. 2nd refinement : Fitting 2D Gaussians to the original image.
- Type:
- neighbor_distance
Mask radius for fitting set to 0.25 of nearest neighbor.
- Type:
0.25
Modules
Atom finding refining
- atomap.atom_finding_refining.calculate_center_of_mass(arr)
Find the center of mass of a NumPy array.
Find the center of mass of a single 2D array, or a list (or multidimensional array) of 2D arrays.
- Parameters:
arr (Numpy 2D Array (or list/nd-array of such))
- Returns:
cy, cx – Giving centre coordinates with sub-pixel accuracy
- Return type:
array of floats (or nd-array of floats)
Examples
>>> import atomap.atom_finding_refining as afr >>> arr = np.random.randint(100, size=(10, 10)) >>> data = afr.calculate_center_of_mass(arr)
Notes
This is a much simpler center of mass approach than the one from scipy. Gotten from stackoverflow: https://stackoverflow.com/questions/37519238/python-find-center-of-object-in-an-image
- atomap.atom_finding_refining.construct_zone_axes_from_sublattice(sublattice, atom_plane_tolerance=0.5, zone_axis_para_list=False, nearest_neighbors=15)
Constructs zone axes for a sublattice.
The zone axes are constructed by finding the 15 nearest neighbors for each atom position in the sublattice, and finding major translation symmetries among the nearest neighbours. Only unique zone axes are kept, and “bad” ones are removed.
- Parameters:
sublattice (Atomap Sublattice)
atom_plane_tolerance (scalar, default 0.5) – When constructing the atomic planes, the method will try to locate the atoms by “jumping” one zone vector, and seeing if there is an atom with the pixel_separation times atom_plane_tolerance. So this value should be increased the atomic planes are non-continuous and “split”.
zone_axis_para_list (parameter list or bool, default False) – A zone axes parameter list is used to name and index the zone axes. See atomap.process_parameters for more info. Useful for automation.
nearest_neighbors (int, default 15) – The number of nearest neighbors which are calculated for each atomic position. Increase this if you want to get more atomic planes.
Example
>>> import atomap.api as am >>> sublattice = am.dummy_data.get_simple_cubic_sublattice() >>> sublattice <Sublattice, (atoms:400,planes:0)> >>> import atomap.atom_finding_refining as afr >>> afr.construct_zone_axes_from_sublattice(sublattice) >>> sublattice <Sublattice, (atoms:400,planes:4)>
Increasing nearest_neighbors to get more planes >>> sublattice = am.dummy_data.get_simple_cubic_sublattice() >>> afr.construct_zone_axes_from_sublattice(sublattice, nearest_neighbors=25) >>> sublattice <Sublattice, (atoms:400,planes:8)>
See also
sublattice._make_translation_symmetry
How unique zone axes are found
sublattice._remove_bad_zone_vectors
How fragmented (“bad zone axis”) are identified and removed.
sublattice.construct_zone_axes
atomap.process_parameters
more info on zone axes parameter list
- atomap.atom_finding_refining.do_pca_on_signal(signal, pca_components=22)
- atomap.atom_finding_refining.find_feature_density(image_data, separation_range=None, separation_step=1, plot_figure=False, plot_debug_figures=False)
Do peak finding with a varying amount of peak separation constrained. Gives a measure of feature density, and what peak separation should be used to find the initial sub-lattice.
Inspiration from the program Smart Align by Lewys Jones.
- atomap.atom_finding_refining.find_features_by_separation(signal, separation_range, separation_step=1, threshold_rel=0.02, pca=False, subtract_background=False, normalize_intensity=False, show_progressbar=True)
Do peak finding with a varying amount of peak separation constrained.
Inspiration from the program Smart Align by Lewys Jones.
- atomap.atom_finding_refining.fit_atom_positions_gaussian(atom_list, image_data, rotation_enabled=True, percent_to_nn=0.4, mask_radius=None, centre_free=True)
Fit a list of Atom_Positions to an image using 2D Gaussians.
The results of the fitting will be saved in the Atom_Position objects themselves, and the old positions will be added to atom.old_pixel_x_list and atom.old_pixel_y_list.
This type of fitting is prone to errors, especially where the atoms are close together, and on noisy images. To reduce the odds of bad fitting the function will rerun the model with a slightly smaller mask_radius or percent_to_nn if: - Any of the fitted positions are outside the model region. - The sigma ratio (highest_sigma/lowest_sigma) is higher than 4.
This repeats 10 times, and if the fitting is still bad, center of mass will be used to set the center position for the atom.
- Parameters:
atom_list (list of Atom_Position objects) – For example the members of a dumbbell.
image_data (NumPy 2D array)
rotation_enabled (bool, optional) – If True (default), the 2D Gaussian will be able to rotate.
percent_to_nn (float, optional)
mask_radius (float, optional) – Radius of the mask around each atom. If this is not set, the radius will be the distance to the nearest atom in the same sublattice times the percent_to_nn value. Note: if mask_radius is not specified, the Atom_Position objects must have a populated nearest_neighbor_list. This is normally done through the sublattice class, but can also be done manually. See below for an example how to do this.
centre_free (bool, optional) – If True (default), the Gaussian will be free to move. Setting this to False can be used to find better values for the other parameters like sigma and A, without the centre positions causing bad fitting.
See also
_make_model_from_atom_list
,_fit_atom_positions_with_gaussian_model
,atom_lattice.Dumbbell_Lattice
,positions
Examples
>>> import numpy as np >>> from atomap.atom_position import Atom_Position >>> from atomap.atom_finding_refining import fit_atom_positions_gaussian
Fitting atomic columns one-by-one
>>> atom_list = [Atom_Position(2, 2), Atom_Position(4, 4)] >>> image = np.zeros((9, 9)) >>> for atom_position in atom_list: ... g_list = fit_atom_positions_gaussian( ... atom_list=[atom_position], image_data=image, mask_radius=2)
Fitting two atoms together
>>> atom_list = [Atom_Position(2, 2), Atom_Position(4, 4)] >>> image = np.zeros((9, 9)) >>> image[2, 2] = 1. >>> image[4, 4] = 1. >>> g_list = fit_atom_positions_gaussian( ... atom_list=atom_list, image_data=image, mask_radius=2)
Not using mask_radius, populating the nearest_neighbor_list manually
>>> image = np.zeros((9, 9)) >>> image[2, 2] = 1. >>> image[5, 5] = 1. >>> atom0 = Atom_Position(2, 2, 0.5, 0.5) >>> atom1 = Atom_Position(5, 5, 0.5, 0.5) >>> atom0.nearest_neighbor_list = [atom1] >>> atom1.nearest_neighbor_list = [atom0] >>> g_list = fit_atom_positions_gaussian([atom0, atom1], image)
- atomap.atom_finding_refining.get_atom_positions(signal, separation=5, threshold_rel=0.02, pca=False, subtract_background=False, normalize_intensity=False, remove_too_close_atoms=True)
Find the most intense features in a HyperSpy signal, where the features has to be separated by a minimum distance.
- Parameters:
signal (HyperSpy 2D signal)
separation (number) – Minimum separation between the features.
threshold_rel (float, default 0.02)
pca (bool, default False) – Do PCA on the signal before doing the peak finding.
subtract_background (bool, default False) – Subtract the average background from the signal before doing the peak finding.
normalize_intensity (bool, default False)
remove_too_close_atoms (bool, default True) – Will attempt to find and remove atoms which are too close to eachother, i.e. less than separation.
- Return type:
NumPy array, list of the most intense atom positions.
Example
>>> import numpy as np >>> import atomap.api as am >>> s = am.dummy_data.get_simple_cubic_signal() >>> from atomap.atom_finding_refining import get_atom_positions >>> atom_positions = get_atom_positions(s, 5) >>> peak_x = atom_positions[:,0] >>> peak_y = atom_positions[:,1]
- atomap.atom_finding_refining.get_feature_separation(signal, separation_range=(5, 30), separation_step=1, pca=False, subtract_background=False, normalize_intensity=False, threshold_rel=0.02, show_progressbar=True)
Plot the peak positions on in a HyperSpy signal, as a function of peak separation.
- Parameters:
signal (HyperSpy signal 2D)
separation_range (tuple, optional, default (5, 30))
separation_step (int, optional, default 1)
pca (bool, default False)
subtract_background (bool, default False)
normalize_intensity (bool, default False)
threshold_rel (float, default 0.02)
show_progressbar (bool, default True)
Example
>>> import numpy as np >>> import hyperspy.api as hs >>> from atomap.atom_finding_refining import get_feature_separation >>> s = hs.signals.Signal2D(np.random.random((500, 500))) >>> s1 = get_feature_separation(s, show_progressbar=False)
- atomap.atom_finding_refining.get_mask_slice(position_list, radius_list)
Get a NumPy slice for extracting atom positions and the region around them.
- atomap.atom_finding_refining.normalize_signal(signal, invert_signal=False)
- atomap.atom_finding_refining.refine_sublattice(sublattice, refinement_config_list, percent_to_nn)
- atomap.atom_finding_refining.subtract_average_background(signal, gaussian_blur=30)
- atomap.atom_finding_refining.zero_array_outside_circle(arr, radius)
Set all values in an array to zero outside a circle defined by radius
Numba jit compatible
- Parameters:
arr (NumPy array) – Must be 2 dimensions
radius (scalar) – Radius of the circle
- Returns:
output_array – Same shape as arr, but with values outside the circle set to zero.
- Return type:
NumPy array
Example
>>> import atomap.atom_finding_refining as afr >>> arr = np.random.randint(100, size=(30, 20))
Analysis tools
- atomap.analysis_tools.get_middle_position_list(sublattice, za0, za1)
Find the middle point between all four neighboring atoms.
The neighbors are found by moving one step along the atom planes belonging to za0 and za1.
So atom planes must be constructed first.
- Parameters:
- Returns:
middle_position_list
- Return type:
Examples
>>> import atomap.analysis_tools as an >>> sublattice = am.dummy_data.get_simple_cubic_sublattice() >>> sublattice.construct_zone_axes() >>> za0 = sublattice.zones_axis_average_distances[0] >>> za1 = sublattice.zones_axis_average_distances[1] >>> middle_position_list = an.get_middle_position_list( ... sublattice, za0, za1)
- atomap.analysis_tools.get_neighbor_middle_position(atom, za0, za1)
Find the middle point between four neighboring atoms.
The neighbors are found by moving one step along the atom planes belonging to za0 and za1.
So atom planes must be constructed first.
- Parameters:
- Returns:
middle_position – If the atom is at the edge by being the last atom in the atom plane, False is returned.
- Return type:
Examples
>>> import atomap.analysis_tools as an >>> sublattice = am.dummy_data.get_simple_cubic_sublattice() >>> sublattice.construct_zone_axes() >>> za0 = sublattice.zones_axis_average_distances[0] >>> za1 = sublattice.zones_axis_average_distances[1] >>> atom = sublattice.atom_list[33] >>> middle_position = an.get_neighbor_middle_position(atom, za0, za1)
- atomap.analysis_tools.get_vector_shift_list(sublattice, position_list)
Find the atom shifts from a central position.
Useful for finding polarization in B-cations in a perovskite structure.
- Parameters:
sublattice (Sublattice object)
position_list (list) – [[x0, y0], [x1, y1], …]
- Returns:
vector_list – In the form [[x0, y0, dx0, dy0], [x1, y1, dx1, dy1]…]
- Return type:
Example
>>> import atomap.analysis_tools as an >>> sublattice = am.dummy_data.get_simple_cubic_sublattice() >>> sublattice.construct_zone_axes() >>> za0 = sublattice.zones_axis_average_distances[0] >>> za1 = sublattice.zones_axis_average_distances[1] >>> middle_position_list = an.get_middle_position_list( ... sublattice, za0, za1) >>> vector_list = an.get_vector_shift_list( ... sublattice, middle_position_list)
- atomap.analysis_tools.pair_distribution_function(image, atom_positions, n_bins=200, rel_range=0.5)
Returns a two dimensional pair distribution function (PDF) from an image of atomic columns.
The intensity of peaks in the PDF is corrected to account for missing information (i.e. the fact atoms are present outside of the field of view) and differences in area at different distances.
- Parameters:
- Returns:
s_pdf – The calculated PDF.
- Return type:
HyperSpy Signal 1D Object
Examples
>>> import atomap.analysis_tools as an >>> s = am.dummy_data.get_simple_cubic_signal() >>> sublattice = am.dummy_data.get_simple_cubic_sublattice() >>> s_pdf = an.pair_distribution_function(s,sublattice.atom_positions) >>> s_pdf.plot()
Quantification
Initial position finding
- class atomap.initial_position_finding.AtomAdderRemover(image, atom_positions=None, distance_threshold=4, norm='linear')
- get_xy_pos_lists()
- is_atom_nearby(x_press, y_press)
- onclick(event)
- replot()
- atomap.initial_position_finding.add_atoms_with_gui(image, atom_positions=None, distance_threshold=4, norm='linear')
Add or remove atoms from a list of atom positions.
Will open a matplotlib figure, where atoms can be added or removed by pressing them.
- Parameters:
image (array-like) – Signal or NumPy array
atom_positions (list of lists, optional) – In the form [[x0, y0], [x1, y1], …]
distance_threshold (int) – Default 4
norm (string, optional) – ‘linear’ or ‘log’, default ‘linear’. If set up log, the intensity will be visualized as a log plot. Useful to see low-intensity atoms.
- Returns:
atom_positions – In the form [[x0, y0], [x1, y1], …]. The list can be updated until the figure is closed.
- Return type:
list of lists
Examples
>>> s = am.dummy_data.get_simple_cubic_signal() >>> atom_positions = am.get_atom_positions(s, separation=9) >>> atom_positions_new = am.add_atoms_with_gui(atom_positions, s)
Log norm
>>> atom_positions_new = am.add_atoms_with_gui( ... atom_positions, s, norm='log')
- atomap.initial_position_finding.find_dumbbell_vector(atom_positions)
Find the vector between the atoms in structure with dumbbells.
For example GaAs-(110).
- Parameters:
atom_positions (list of lists) – In the form [[x0, y0], [x1, y1], …]
- Returns:
dumbbell_vector
- Return type:
Examples
>>> import atomap.initial_position_finding as ipf >>> s = am.dummy_data.get_dumbbell_signal() >>> atom_positions = am.get_atom_positions(s, separation=4) >>> dumbbell_vector = ipf.find_dumbbell_vector(atom_positions)
- atomap.initial_position_finding.make_atom_lattice_dumbbell_structure(s, dumbbell_positions, dumbbell_vector, show_progressbar=True)
Make Atom_Lattice object from image of dumbbell structure.
- Parameters:
- Returns:
dumbbell_lattice
- Return type:
Atomap Dumbbell_Lattice object
Examples
>>> import atomap.initial_position_finding as ipf >>> s = am.dummy_data.get_dumbbell_signal() >>> dumbbell_positions = am.get_atom_positions(s, separation=16) >>> atom_positions = am.get_atom_positions(s, separation=4) >>> dumbbell_vector = ipf.find_dumbbell_vector(atom_positions) >>> dumbbell_lattice = ipf.make_atom_lattice_dumbbell_structure( ... s, dumbbell_positions, dumbbell_vector)
- atomap.initial_position_finding.select_atoms_with_gui(image, atom_positions, verts=None, invert_selection=False)
Get a subset of a list of atom positions.
Will open a matplotlib figure, where a polygon is drawn interactively.
- Parameters:
image (array-like) – Signal or NumPy array
atom_positions (list of lists, or NumPy array) – In the form [[x0, y0], [x1, y1], …]
verts (list of tuples) – List of positions, spanning an enclosed region. [(x0, y0), (x1, y1), …]. Need to have at least 3 positions.
invert_selection (bool, optional) – Get the atom positions outside the region, instead of the ones inside it. Default False.
- Returns:
atom_positions_selected
- Return type:
list of lists
Examples
>>> s = am.dummy_data.get_simple_cubic_signal() >>> atom_positions = am.get_atom_positions(s, separation=9) >>> atom_positions_selected = am.select_atoms_with_gui( ... s, atom_positions) >>> sublattice = am.Sublattice(atom_positions_selected, s)
Using the function non-interactively
>>> verts = [(60, 235), (252, 236), (147, 58)] >>> atom_positions_selected = am.select_atoms_with_gui( ... s, atom_positions, verts=verts) >>> sublattice = am.Sublattice(atom_positions_selected, s)
Getting the atom positions outside the selected regions, instead of outside it
>>> atom_positions_selected = am.select_atoms_with_gui( ... s, atom_positions, verts=verts, invert_selection=True) >>> sublattice = am.Sublattice(atom_positions_selected, s)
Input output (IO)
Module containing functions to save and load Atom_Lattice objects.
- atomap.io.load_atom_lattice_from_hdf5(filename, construct_zone_axes=True)
Load an Atomap HDF5-file, restoring a saved Atom_Lattice.
- Parameters:
filename (string) – Filename of the HDF5-file.
construct_zone_axes (bool) – If True, find relations between atomic positions by constructing atomic planes. Default True.
- Return type:
Atomap Atom_Lattice object
Main
- atomap.main.make_atom_lattice_from_image(s_image0, process_parameter=None, pixel_separation=None, s_image1=None, debug_plot=False)
- atomap.main.run_image_filtering(signal, invert_signal=False)
Subtracts background, filters noise with PCA, and normalizes a signal.
- Parameters:
signal (HyperSpy Signal2D)
invert_signal (bool, default False) – Inverts the image to 1./signal.data
- Returns:
filtered_signal
- Return type:
HyperSpy Signal2D
Example
>>> import atomap.api as am >>> s = am.dummy_data.get_simple_cubic_signal() >>> from atomap.main import run_image_filtering >>> s_new = run_image_filtering(s)
Plotting
- atomap.plotting.get_rgb_array(angle, magnitude, rotation=0, angle_lim=None, magnitude_lim=None)
- atomap.plotting.make_color_wheel(ax, rotation=0)
- atomap.plotting.normalize_array(np_array, max_number=1.0)
- atomap.plotting.plot_complex_image_map_line_profile_using_interface_plane(image_data, amplitude_data, phase_data, line_profile_amplitude_data, line_profile_phase_data, interface_plane, atom_plane_list=None, data_scale=1, atom_list=None, extra_marker_list=None, amplitude_image_lim=None, phase_image_lim=None, plot_title='', add_color_wheel=False, color_bar_markers=None, vector_to_plot=None, rotate_atom_plane_list_90_degrees=False, prune_outer_values=False, figname='map_data.jpg')
- Parameters:
atom_list (list of Atom_Position instances)
extra_marker_list (two arrays of x and y [[x_values], [y_values]])
- atomap.plotting.plot_feature_density(separation_value_list, peakN_list, figname='feature_density.png')
- atomap.plotting.plot_image_map_line_profile_using_interface_plane(image_data, heatmap_data_list, line_profile_data_list, interface_plane, atom_plane_list=None, data_scale=1, data_scale_z=1, atom_list=None, extra_marker_list=None, clim=None, plot_title='', vector_to_plot=None, rotate_atom_plane_list_90_degrees=False, prune_outer_values=False, figname='map_data.jpg')
- Parameters:
atom_list (list of Atom_Position instances)
extra_marker_list (two arrays of x and y [[x_values], [y_values]])
- atomap.plotting.plot_line_profiles_from_parameter_input(parameter_list, parameter_name_list=None, invert_line_profiles=False, extra_line_marker_list=None, x_lim=False, figname='line_profile_list.jpg')
- atomap.plotting.plot_vector_field(x_pos_list, y_pos_list, x_rot_list, y_rot_list, save=True)
- atomap.plotting.plot_zone_vector_and_atom_distance_map(image_data, distance_data, atom_planes=None, distance_data_scale=1, atom_list=None, extra_marker_list=None, clim=None, atom_plane_marker=None, plot_title='', vector_to_plot=None, figsize=(10, 20), figname='map_data.jpg')
- Parameters:
atom_list (list of Atom_Position instances)
extra_marker_list (two arrays of x and y [[x_values], [y_values]])
- atomap.plotting.vector_list_to_marker_list(vector_list, color='red', scale=1.0)
Make a marker list from a vector list.
- Parameters:
vector_list (list) – In the form [[x0, y0, dx0, dy0], …]
color (string, optional) – Default ‘red’
scale (scalar, optional) – Default 1.
- Returns:
marker_list
- Return type:
list of markers
Examples
>>> import atomap.plotting as pl >>> vector_list = [[13, 11, -2, 1], [20, 12, 2, -3]] >>> marker_list = pl.vector_list_to_marker_list( ... vector_list, color='red', scale=1.)
Tools
- class atomap.tools.Fingerprinter(cluster_algo=DBSCAN(eps=0.1, min_samples=10))
Produces a distance-fingerprint from an array of neighbor distance vectors.
To avoid introducing our own interface we’re going to use scikit-learn Estimator conventions to name the method, which produces our fingerprint, ‘fit’ and store our estimations as attributes with a trailing underscore in their names. http://scikit-learn.org/stable/developers/contributing.html#fitting http://scikit-learn.org/stable/developers/contributing.html#estimated-attributes)
- fingerprint_
The main result. The contents of fingerprint_ can be described as the relative distances to neighbours in the generalized neighborhood.
- Type:
array, shape = (n_clusters,)
- cluster_centers_
The cluster center coordinates from which the fingerprint was produced.
- Type:
array, shape = (n_clusters, n_dimensions)
- cluster_algo_.labels_
Integer labels that denote which cluster each point belongs to.
- Type:
array, shape = (n_points,)
- fit(X, max_neighbors=150000)
- Parameters:
X (array, shape = (n_points, n_dimensions)) – This array is typically a transpose of a subset of the returned value of sublattice.get_nearest_neighbor_directions_all()
max_neighbors (int, default 150000) – If the length of X is larger than max_neighbors, X will be reduced to max_neighbors. The selection is random. This is done to allow very large datasets to be processed, since having X too large causes the fitting to use too much memory.
Notes
More information about memory use: http://scikit-learn.org/stable/modules/clustering.html#dbscan
- atomap.tools.array2signal1d(array, scale=1.0, offset=0.0)
- atomap.tools.array2signal2d(numpy_array, scale=1.0, units='', rotate_flip=False)
- atomap.tools.calculate_angle(v1, v2)
- atomap.tools.calculate_point_record(image, points, max_radius)
Creates a Voronoi array where equal values belong to the same Voronoi cell
- Parameters:
point_record (2D zero array of same shape as the image to be mapped)
points (Array like of shape (2,N))
max_radius (Integer, max radius of each Voronoi Cell)
- Returns:
point_record (Voronoi array where equal values belong to)
the same Voronoi cell
- atomap.tools.combine_projected_positions_layers(layer_list)
Get the mean position of a projected property, and its position.
Will also include the standard deviation in of the positions in the layers.
- atomap.tools.dotproduct(v1, v2)
- atomap.tools.find_atom_position_between_atom_planes(image, atom_plane0, atom_plane1, orthogonal_zone_vector, integration_width_percent=0.2, max_oxygen_sigma_percent=0.2)
- atomap.tools.find_atom_positions_for_all_atom_planes(image, sublattice, parallel_zone_vector, orthogonal_zone_vector)
- atomap.tools.find_atom_positions_for_an_atom_plane(image, atom_plane0, atom_plane1, orthogonal_zone_vector)
- atomap.tools.find_average_distance_between_atoms(input_data_list, crop_start=3, crop_end=3, threshold=0.4)
Return the distance between monolayers.
Returns the maximal separation between two adjacent points in input_data_list[:, 0], as a good approximation for monolayer separation.
- Parameters:
input_data_list (NumPy array) – An array where the distance from a point to a line is in input_data_list[:, 0]
crop_start (int) – Before and after the index given by crop_start and crop_end, the data is ignored. By default 3, to ignore outliers.
crop_end (int) – Before and after the index given by crop_start and crop_end, the data is ignored. By default 3, to ignore outliers.
threshold (float, default 0.4) – Atoms with a separation of more than the threshold times the largest atom separation will be counted as members of different planes.
- Returns:
first_peak (float) – The monolayer separation.
monolayer_sep (array) – An array with monolayer separations
mean_separation (float) – The mean monolayer separation
Example
>>> import atomap.tools as to >>> position_list = [] >>> for i in range(0, 100, 9): ... positions = np.ones(10) * i ... position_list.extend(positions) >>> position_list = np.array(position_list) >>> property_list = np.ones(len(position_list)) >>> input_data_list = np.stack((position_list, property_list), axis=1) >>> output = to.find_average_distance_between_atoms(input_data_list) >>> first_peak, monolayer_sep, mean_separation = output
- atomap.tools.find_smallest_distance(i, j, points)
Finds the smallest distance between coordinates (i, j) and a list of coordinates.
- Parameters:
i (Integer)
j (Integer)
points (array like of shape (2,N))
- Returns:
distMin (Minimum distance)
minIndex (Index of minimum distance in points)
Example
>>> import numpy as np >>> points = np.random.random((2, 10000)) >>> i, j = 0.5, 0.5 >>> smallest_distance = find_smallest_distance(i, j, points)
- atomap.tools.fliplr_points_and_signal(signal, x_array, y_array)
Horizontally flip a set of points and a HyperSpy signal.
For making sure both the image and points are flipped correctly.
- Parameters:
signal (HyperSpy 2D signal)
x_array (array-like)
y_array (array-like)
- Returns:
flipped_signal (HyperSpy signal)
flipped_x_array, flipped_y_array (NumPy arrays)
Examples
>>> import atomap.api as am >>> import atomap.tools as to >>> sublattice = am.dummy_data.get_distorted_cubic_sublattice() >>> s = sublattice.get_atom_list_on_image() >>> x, y = sublattice.x_position, sublattice.y_position >>> s_flip, x_flip, y_flip = to.fliplr_points_and_signal(s, x, y)
Plotting the data
>>> import matplotlib.pyplot as plt >>> fig, ax = plt.subplots() >>> cax = ax.imshow(s_flip.data, origin='lower', ... extent=s_flip.axes_manager.signal_extent) >>> cpoints = ax.scatter(x_flip, y_flip)
- atomap.tools.fliplr_points_around_signal_centre(signal, x_array, y_array)
Horizontally flip a set of points around the centre of a signal.
- Parameters:
signal (HyperSpy 2D signal)
x_array (array-like)
y_array (array-like)
- Returns:
flipped_x_array, flipped_y_array
- Return type:
NumPy arrays
Examples
>>> import atomap.api as am >>> import atomap.tools as to >>> sublattice = am.dummy_data.get_distorted_cubic_sublattice() >>> s = sublattice.get_atom_list_on_image() >>> x, y = sublattice.x_position, sublattice.y_position >>> x_rot, y_rot = to.fliplr_points_around_signal_centre(s, x, y)
- atomap.tools.get_arbitrary_slice(image, start_point, end_point, width, debug_figname=None)
- atomap.tools.get_atom_planes_square(sublattice, atom_plane1, atom_plane2, interface_atom_plane, zone_vector, debug_plot=False)
- Parameters:
sublattice (Atomap sublattice object)
atom_plane1 (Atomap atom_plane object)
atom_plane2 (Atomap atom_plane object)
- atomap.tools.get_integrated_intensity(point_record, image, point_index, include_edge_cells=True)
Using a Voronoi point_record array, integrate a (minimum 2D) image array at each pixel
- Parameters:
point_record (2D zero array of same shape as the image to be mapped)
image (The ndarray to integrate the voronoi cells on)
point_index (Array like of shape (2,N))
- Returns:
integrated_record (Voronoi array where equal values belong to)
the same Voronoi cell
- atomap.tools.get_point_between_four_atoms(atom_list)
Get the mean point between four atom position objects.
- Parameters:
atom_list (list) – List of four atom position objects
- Returns:
middle_position – (x_pos, y_pos)
- Return type:
Example
>>> import atomap.tools as at >>> import atomap.atom_position as ap >>> atom0, atom1 = ap.Atom_Position(10, 30), ap.Atom_Position(20, 30) >>> atom2, atom3 = ap.Atom_Position(10, 40), ap.Atom_Position(20, 40) >>> mid_pos = at.get_point_between_four_atoms((atom0, atom1, atom2, atom3))
- atomap.tools.get_point_between_two_atoms(atom_list)
- atomap.tools.get_slice_between_four_atoms(image, start_atoms, end_atoms, width)
- atomap.tools.get_slice_between_two_atoms(image, atom0, atom1, width)
- atomap.tools.integrate(s, points_x, points_y, method='Voronoi', max_radius='Auto', show_progressbar=True, remove_edge_cells=False, edge_pixels=1)
Given a spectrum image a set of points and a maximum outer radius, this function integrates around each point in an image, using either Voronoi cell or watershed segmentation methods.
- Parameters:
s (HyperSpy signal or array-like object) – Assuming 2D, 3D or 4D dataset where the spatial dimensions are 2D and any remaining dimensions are spectral.
point_x (list) – Detailed list of the x and y coordinates of each point of interest within the image.
point_y (list) – Detailed list of the x and y coordinates of each point of interest within the image.
method (string) – ‘Voronoi’ or ‘Watershed’
max_radius ({'Auto'} int) – A maximum outer radius for each Voronoi Cell. If a pixel exceeds this radius it will not be included in the cell. This allows analysis of a surface and particles. If ‘max_radius’ is left as ‘Auto’ then it will be set to the largest dimension in the image.
remove_edge_cells (bool) – Determine whether to replace the cells touching the signal edge with np.nan values, which makes automatic contrast estimation easier later
edge_pixels (int) – Only used if remove_edge_cells is True. Determines the number of pixels from the border to remove.
show_progressbar (bool, optional) – Default True
- Returns:
integrated_intensity (NumPy array) – An array where dimension 0 is the same length as points, and subsequent subsequent dimension are energy dimensions.
intensity_record (HyperSpy signal, same size as s) – Each pixel/voxel in a particular segment or region has the value of the integration, value.
point_record (HyperSpy Signal2D with no navigation dimension) – Image showing where each integration region is, pixels equating to point 0 (integrated_intensity[0]) all have value 0, all pixels equating to integrated_intensity[1] all have value 1 etc.
Examples
>>> import atomap.api as am >>> from atomap.tools import integrate >>> import hyperspy.api as hs >>> sublattice = am.dummy_data.get_simple_cubic_sublattice( ... image_noise=True) >>> image = hs.signals.Signal2D(sublattice.image) >>> i_points, i_record, p_record = integrate( ... image, ... points_x=sublattice.x_position, ... points_y=sublattice.y_position, method='Voronoi') >>> i_record.plot()
For a 3 dimensional dataset, with artificial EELS data
>>> s = am.dummy_data.get_eels_spectrum_survey_image() >>> s_eels = am.dummy_data.get_eels_spectrum_map() >>> peaks = am.get_atom_positions(s, separation=4) >>> i_points, i_record, p_record = integrate( ... s_eels, peaks[:, 0], peaks[:, 1], max_radius=3)
Note
Works in principle with 3D and 4D data sets but will quickly hit a memory error with large sizes.
- atomap.tools.length(v)
- atomap.tools.project_position_property(input_data_list, interface_plane)
Project 2D positions onto a 1D plane.
The 2D positions are found as function of distance to the interface_plane. In this case, one will get the positions as a function of atomic plane from the interface_plane.
Note, positions will not be projected on to the interface_plane, but on a plane perpendicular to the interface_plane. The returned data will give the distance from the projected positions to the interface_plane.
- Parameters:
input_data_list (Numpy array, [Nx3]) – Numpy array with positions and property value. Must be in the from [[x,y,z]], so that the x-positions are extracted using input_data_list[:,0]. y-positions using input_data_list[:,1]. Property value using input_data_list[:,2].
interface_plane (Atomap atom_plane object)
- Returns:
Data list – Array in the form [[position, property]].
- Return type:
NumPy Array
Example
>>> from numpy.random import random >>> from atomap.sublattice import Sublattice >>> pos = [[x, y] for x in range(9) for y in range(9)] >>> sublattice = Sublattice(pos, random((9, 9))) >>> sublattice.construct_zone_axes() >>> x, y = sublattice.x_position, sublattice.y_position >>> z = sublattice.ellipticity >>> input_data_list = np.array([x, y, z]).swapaxes(0, 1) >>> from atomap.tools import project_position_property >>> plane = sublattice.atom_plane_list[10] >>> data = project_position_property(input_data_list, plane) >>> positions = data[:,0] >>> property_values = data[:,1] >>> cax = plt.plot(positions, property_values)
- atomap.tools.remove_atoms_from_image_using_2d_gaussian(image, sublattice, percent_to_nn=0.4, show_progressbar=True)
- Parameters:
image (NumPy 2D array)
sublattice (Atomap sublattice object)
percent_to_nn (float) – Percent to nearest neighbor. The function will find the closest nearest neighbor to the current atom position, and this value times percent_to_nn will be the radius of the mask centered on the atom position. Value should be somewhere between 0.01 (1%) and 1 (100%). Having a too low value might lead to bad fitting.
show_progressbar (bool, default True)
- Returns:
subtracted_image
- Return type:
NumPy 2D array
Examples
>>> atom_lattice = am.dummy_data.get_simple_atom_lattice_two_sublattices() >>> sublattice0 = atom_lattice.sublattice_list[0] >>> sublattice0.find_nearest_neighbors() >>> import atomap.tools as at >>> image_subtracted = at.remove_atoms_from_image_using_2d_gaussian( ... image=atom_lattice.image, sublattice=sublattice0, ... show_progressbar=False) >>> import hyperspy.api as hs >>> s = hs.signals.Signal2D(image_subtracted) >>> s.plot()
Decrease percent_to_nn, to reduce the effect of overlapping atoms. For this dataset it won’t change much, but might be very useful for real datasets.
>>> image_subtracted = at.remove_atoms_from_image_using_2d_gaussian( ... image=atom_lattice.image, sublattice=sublattice0, ... percent_to_nn=0.2, show_progressbar=False)
- atomap.tools.remove_integrated_edge_cells(i_points, i_record, p_record, edge_pixels=1, use_nans=True, inplace=False)
Removes any cells that touch within a number of pixels of the image border.
Note on using use_nans: If this is used on a dataset with more than two dimensions, the resulting HyperSpy i_record signal might be needed to be viewed with i_record.plot(navigator=’slider’), since HyperSpy may throw an error when plotting a dataset with only NaNs present.
- Parameters:
i_points (NumPy array) – The output of the Atomap integrate function or method
i_record (HyperSpy signal) – The output of the Atomap integrate function or method
p_record (HyperSpy signal) – The output of the Atomap integrate function or method
- Returns:
i_points (NumPy array) – Modified list of integrated intensities with either np.nan or 0 on the removed values, which preserves the atom index.
i_record (HyperSpy signal) – Modified integrated intensity record, with either np.nan or 0 on the removed values, which preserves the atom index
p_record (HyperSpy signal, same size as image) – Modified points record, where removed areas have value = -1.
Example
>>> s = am.dummy_data.get_fantasite() >>> points_x, points_y = am.get_atom_positions(s).T >>> i, ir, pr = am.integrate( ... s, ... points_x, ... points_y, ... method='Voronoi', ... remove_edge_cells=False) >>> from atomap.tools import remove_integrated_edge_cells >>> i2, ir2, pr2 = remove_integrated_edge_cells( ... i, ir, pr, edge_pixels=5, use_nans=True)
- atomap.tools.rotate_points_and_signal(signal, x_array, y_array, rotation)
Rotate a set of points and a HyperSpy signal.
For making sure both the image and points are rotated correctly.
- Parameters:
signal (HyperSpy 2D signal)
x_array (array-like)
y_array (array-like)
rotation (scalar) – In degrees
- Returns:
rotated_signal (HyperSpy signal)
rotated_x_array, rotated_y_array (NumPy arrays)
Examples
>>> import atomap.api as am >>> import atomap.tools as to >>> sublattice = am.dummy_data.get_distorted_cubic_sublattice() >>> s = sublattice.get_atom_list_on_image() >>> x, y = sublattice.x_position, sublattice.y_position >>> s_rot, x_rot, y_rot = to.rotate_points_and_signal(s, x, y, 30)
Plotting the data
>>> fig, ax = plt.subplots() >>> cax = ax.imshow(s_rot.data, origin='lower', ... extent=s_rot.axes_manager.signal_extent) >>> cpoints = ax.scatter(x_rot, y_rot)
- atomap.tools.rotate_points_around_signal_centre(signal, x_array, y_array, rotation)
Rotate a set of points around the centre of a signal.
- Parameters:
signal (HyperSpy 2D signal)
x_array (array-like)
y_array (array-like)
rotation (scalar) – In degrees
- Returns:
rotated_x_array, rotated_y_array
- Return type:
NumPy arrays
Examples
>>> import atomap.api as am >>> import atomap.tools as to >>> sublattice = am.dummy_data.get_distorted_cubic_sublattice() >>> s = sublattice.get_atom_list_on_image() >>> x, y = sublattice.x_position, sublattice.y_position >>> x_rot, y_rot = to.rotate_points_around_signal_centre(s, x, y, 30)
- atomap.tools.sort_positions_into_layers(data_list, layer_distance, combine_layers=True)
Combine clustered positions into groups.
Atoms with a similar distance from a line belong to the same plane parallel to this line. Atoms in data_list are grouped based on which plane they belong to. If there is only one atom in a layer, it will be disregarded as it gives a high uncertainty.
- Parameters:
data_list (NumPy array) – An array where the distance from a point to a line is in input_data_list[:, 0], and the property of the point (atom) is in [:,1]
layer_distance (float) – The half width of a layer, used to determine which layer an atom belongs to.
combine_layers (bool, default True) – If True, the values for distance and property is averaged for each layer.
- Returns:
layer_list – If combine_layers is True, a list of the average position and property of the points in the layer. If False, a nested list where each element in layer_list contains a list the atoms (position and property) in the layer.
- Return type:
- atomap.tools.sort_projected_positions_into_layers(projected_positions, margin=0.5)
Dummy data
- atomap.dummy_data.get_atom_counting_signal()
Get an artificial image of 4 different atomic column intensities.
- Returns:
s_atom_counting
- Return type:
HyperSpy 2D signal
Examples
>>> import atomap.api as am >>> s = am.dummy_data.get_atom_counting_signal() >>> s.plot()
- atomap.dummy_data.get_distorted_cubic_signal(image_noise=False)
Generate a test image signal of a distorted cubic atomic structure.
- Parameters:
image_noise (default False) – If True, will add Gaussian noise to the image.
- Returns:
signal
- Return type:
HyperSpy 2D
Examples
>>> import atomap.api as am >>> s = am.dummy_data.get_distorted_cubic_signal() >>> s.plot()
- atomap.dummy_data.get_distorted_cubic_sublattice(image_noise=False)
Generate a test sublattice of a distorted cubic atomic structure.
- Parameters:
image_noise (default False) – If True, will add Gaussian noise to the image.
- Returns:
sublattice
- Return type:
Atomap Sublattice
Examples
>>> import atomap.api as am >>> sublattice = am.dummy_data.get_distorted_cubic_sublattice() >>> sublattice.plot()
- atomap.dummy_data.get_dumbbell_heterostructure_dumbbell_lattice()
Get a dumbbell heterostructure dumbbell lattice.
- Returns:
dumbbell_lattice
- Return type:
Atomap Dumbbell_Lattice
Example
>>> dl = am.dummy_data.get_dumbbell_heterostructure_dumbbell_lattice() >>> dl.plot()
- atomap.dummy_data.get_dumbbell_heterostructure_signal()
Get a dumbbell heterostructure signal.
- Returns:
s_dumbbell
- Return type:
HyperSpy Signal2D
Example
>>> s = am.dummy_data.get_dumbbell_heterostructure_signal() >>> s.plot()
- atomap.dummy_data.get_dumbbell_signal()
Get a simple dumbbell signal.
- Returns:
s_dumbbell
- Return type:
HyperSpy Signal2D
Example
>>> s = am.dummy_data.get_dumbbell_signal() >>> s.plot()
- atomap.dummy_data.get_eels_spectrum_map(add_noise=True)
Get an artificial atomic resolution EELS map of LaMnO3
Containing the Mn-L23 and La-M54 edges.
- Parameters:
add_noise (bool) – If True, will add Gaussian noise to the spectra. Default True.
- Returns:
eels_map
- Return type:
HyperSpy EELSSpectrum
Example
>>> import atomap.api as am >>> s_eels_map = am.dummy_data.get_eels_spectrum_map() >>> s_eels_map.plot()
Not adding noise
>>> s_eels_map = am.dummy_data.get_eels_spectrum_map(add_noise=False)
See also
get_eels_spectrum_survey_image
signal with same spatial dimensions
- atomap.dummy_data.get_eels_spectrum_survey_image()
Get an artificial survey image of atomic resolution EELS map of LaMnO3
- Returns:
survey_image
- Return type:
HyperSpy Signal2D
Example
>>> import atomap.api as am >>> s = am.dummy_data.get_eels_spectrum_survey_image() >>> s.plot()
See also
get_eels_spectrum_map
corresponding EELS map
- atomap.dummy_data.get_fantasite()
Fantasite is a fantastic structure with several interesting structural variations.
It contains two sublattices, domains with elliptical atomic columns and tilt-patterns. This function returns a HyperSpy 2D signal.
Examples
>>> import atomap.api as am >>> s = am.dummy_data.get_fantasite() >>> s.plot()
- Returns:
fantasite_signal
- Return type:
HyperSpy 2D signal
See also
get_fantasite_sublattice
get a sublattice object of the fantasite.
- atomap.dummy_data.get_fantasite_atom_lattice()
Fantasite is a fantastic structure with several interesting structural variations.
It contains two sublattices, domains with elliptical atomic columns and tilt-patterns. This function returns an Atomap Atom_Lattice.
Examples
>>> import atomap.api as am >>> atom_lattice = am.dummy_data.get_fantasite_atom_lattice() >>> atom_lattice.plot()
- atomap.dummy_data.get_fantasite_sublattice()
Fantasite is a fantastic structure with several interesting structural variations.
It contains two sublattices, domains with elliptical atomic columns and tilt-patterns. This function returns an Atomap sublattice.
Currently this function returns the two sublattices as one sublattice. To get these sublattices separately, see get_fantasite_atom_lattice
Examples
>>> import atomap.api as am >>> sublattice = am.dummy_data.get_fantasite_sublattice() >>> sublattice.plot()
See also
get_fantasite
get a Signal2D object of the fantasite.
get_fantasite_atom_lattice
get the atom lattice for fantasite, with both sublattices.
- atomap.dummy_data.get_hexagonal_double_signal(image_noise=False)
Generate a test image signal of a hexagonal structure.
Similar to MoS2.
- Parameters:
image_noise (default False) – If True, will add Gaussian noise to the image.
- Returns:
>>> s = am.dummy_data.get_hexagonal_double_signal()
>>> s.plot()
Adding image noise
>>> s = am.dummy_data.get_hexagonal_double_signal(image_noise=True)
>>> s.plot()
- atomap.dummy_data.get_nanoparticle_signal()
Get an simulated nanoparticle image.
- Returns:
signal
- Return type:
HyperSpy Signal2D
Example
>>> import atomap.api as am >>> s = am.dummy_data.get_nanoparticle_signal() >>> s.plot()
- atomap.dummy_data.get_perovskite110_ABF_signal(image_noise=False)
Returns signals representing an ABF image of a perovskite <110>.
- Parameters:
image_noise (default False) – If True, will add Gaussian noise to the image.
- Returns:
signal
- Return type:
HyperSpy 2D,
Examples
>>> import atomap.api as am >>> s_ABF = am.dummy_data.get_perovskite110_ABF_signal()
- atomap.dummy_data.get_perovskite_001_atom_lattice(image_noise=False, set_element_info=True)
- atomap.dummy_data.get_perovskite_001_signal(image_noise=False)
- atomap.dummy_data.get_polarization_film_atom_lattice(image_noise=False)
Get an Atom Lattice with two sublattices resembling a perovskite film.
Similar to a perovskite oxide thin film, where the B cations are shifted in the film.
- Parameters:
image_noise (bool, default False)
- Returns:
simple_atom_lattice
- Return type:
Atom_Lattice object
Examples
>>> import atomap.api as am >>> atom_lattice = am.dummy_data.get_polarization_film_atom_lattice() >>> atom_lattice.plot()
With image noise
>>> atom_lattice = am.dummy_data.get_polarization_film_atom_lattice( ... image_noise=True) >>> atom_lattice.plot()
- atomap.dummy_data.get_polarization_film_signal(image_noise=False)
Get a signal with two sublattices resembling a perovskite film.
Similar to a perovskite oxide thin film, where the B cations are shifted in the film.
- Parameters:
image_noise (bool, optional) – Default False
- Returns:
signal
- Return type:
HyperSpy 2D signal
Examples
>>> import atomap.api as am >>> s = am.dummy_data.get_polarization_film_signal() >>> s.plot()
With image noise
>>> s = am.dummy_data.get_polarization_film_signal(image_noise=True) >>> s.plot()
- atomap.dummy_data.get_precipitate_signal()
Get a signal emulating a preciptate embedded in a matrix.
- Returns:
signal_percipitate
- Return type:
HyperSpy Signal2D
Examples
>>> s = am.dummy_data.get_precipitate_signal() >>> s.plot()
- atomap.dummy_data.get_scanning_distortion_sublattice()
Get an artificial sublattice emulating scanning distortions.
- Returns:
sublattice
- Return type:
Atomap Sublattice
Example
>>> import atomap.api as am >>> sublattice = am.dummy_data.get_scanning_distortion_sublattice() >>> sublattice.plot()
- atomap.dummy_data.get_simple_atom_lattice_two_sublattices(image_noise=False)
Returns a simple atom_lattice with two sublattices.
- Parameters:
image_noise (bool, default False)
- Returns:
simple_atom_lattice
- Return type:
Atom_Lattice object
Examples
>>> import atomap.api as am >>> al = am.dummy_data.get_simple_atom_lattice_two_sublattices() >>> al.plot()
- atomap.dummy_data.get_simple_cubic_signal(image_noise=False)
Generate a test image signal of a simple cubic atomic structure.
- Parameters:
image_noise (default False) – If True, will add Gaussian noise to the image.
- Returns:
signal
- Return type:
HyperSpy 2D
Examples
>>> import atomap.api as am >>> s = am.dummy_data.get_simple_cubic_signal() >>> s.plot()
- atomap.dummy_data.get_simple_cubic_sublattice(image_noise=False)
Generate a test sublattice of a simple cubic atomic structure.
- Parameters:
image_noise (default False) – If True, will add Gaussian noise to the image.
- Returns:
sublattice
- Return type:
Atomap Sublattice
Examples
>>> import atomap.api as am >>> sublattice = am.dummy_data.get_simple_cubic_sublattice() >>> sublattice.plot()
- atomap.dummy_data.get_simple_cubic_with_vacancies_signal(image_noise=False)
Generate a test image signal with some vacancies
Same as the simple cubic signal, but with 3 missing atoms.
- Parameters:
image_noise (default False) – If True, will add Gaussian noise to the image.
- Returns:
signal
- Return type:
HyperSpy 2D
Examples
>>> import atomap.api as am >>> s = am.dummy_data.get_simple_cubic_with_vacancies_signal() >>> s.plot()
With image noise
>>> s1 = am.dummy_data.get_simple_cubic_with_vacancies_signal( ... image_noise=True) >>> s1.plot()
- atomap.dummy_data.get_simple_cubic_with_vacancies_sublattice(image_noise=False)
Generate a test sublattice with some vacancies
Same as the simple cubic sublattice, but with 3 missing atoms.
- Parameters:
image_noise (default False) – If True, will add Gaussian noise to the image.
- Returns:
signal
- Return type:
HyperSpy 2D
Examples
>>> import atomap.api as am >>> s = am.dummy_data.get_simple_cubic_with_vacancies_sublattice() >>> s.plot()
With image noise
>>> s1 = am.dummy_data.get_simple_cubic_with_vacancies_sublattice( ... image_noise=True) >>> s1.plot()
- atomap.dummy_data.get_simple_heterostructure_signal(image_noise=True)
- atomap.dummy_data.get_single_atom_sublattice()
Get a sublattice containing a single atom.
- Returns:
sublattice_single_atom
- Return type:
Atomap Sublattice class
Example
>>> sublattice = am.dummy_data.get_single_atom_sublattice() >>> sublattice.plot()
- atomap.dummy_data.get_two_sublattice_signal()