The APyT spatial distribution map (SDM) module

Introduction

Atom probe tomography (APT) reconstructions are known to introduce geometric distortions, particularly affecting the fidelity of the crystal lattice. While the resolution along the detection (\(z\)) axis is typically sufficient to resolve individual lattice planes—and to determine the correct \(z\)-scaling—the scaling of the lateral (\(x\) and \(y\)) directions is less straightforward.

Geiser et al. demonstrated that valuable crystallographic information can still be extracted from lateral directions using spatial distribution maps (SDMs).

An SDM is a two-dimensional histogram of interatomic distance vectors \(\Delta \vec{r}_{ij} = (\Delta x_{ij}, \Delta y_{ij}, \Delta z_{ij})^T\), where the lateral components \(\Delta x_{ij}\) and \(\Delta y_{ij}\) are used as histogram axes. Depending on the crystal structure and orientation, certain combinations \((\Delta x_{ij}, \Delta y_{ij})\) will occur more frequently, appearing as distinct maxima in the histogram.

General Procedure

As shown by Geiser et al., the accuracy of lateral spatial resolution improves when atomic pairs are selected within a narrow separation along the \(z\)-axis. This requires optimal alignment of the sample such that the \(z\)-axis corresponds closely to a crystallographic direction. When constructed under these conditions, SDMs reveal distinct (but distorted) maxima.

Since the ideal crystallographic positions of these maxima are known, a system of linear equations can be formulated and solved to transform the distorted SDM into its correct geometric representation. This transformation also allows for the determination of the remaining \(z\)-scaling factor.

List of functions

This module provides several functions for generating, analyzing, and rectifying SDMs from 3D APT data. These include tools for alignment optimization, SDM construction, and spatial correction.

  • lattice_vectors(): Get vectors which span the lattice planes.

  • optimize_alignment(): Get optimal alignment in normal direction.

  • rdf(): Analytic model of (one-dimensional) radial distribution function in normal direction.

  • rdf_1d(): Create histogram of radial distribution function in normal direction.

  • rdf_lateral(): Create histogram of radial distribution function for lateral directions.

  • rectify_sdm(): Rectify SDMs (and atomic positions).

  • rotate(): Rotate positions around two axes for alignment in normal direction.

  • sdm(): Create SDMs.

apyt.analysis.crystallography.sdm.lattice_vectors(data, angles, dir_key)[source]

Get vectors which span the lattice planes before and after rotation.

If the lattice planes are not perfectly aligned in normal direction, the box dimensions in lateral directions do not reflect the actual cross-sectional area covered by the lattice planes due to the inclination.

This method accounts for this specific inclination and returns the vectors which span the lattice planes in lateral directions before and after the angular alignment (see optimize_alignment()). These vectors can then be used to calculate the actual cross-sectional area of the lattice planes with the aid of the cross product.

Note that the value obtained by this correction may not differ significantly from the lateral box dimensions if the alignment of the lattice planes in normal direction is already sufficiently satisfied (i.e. very small Euler angles).

Parameters:
  • data (ndarray, shape (n, 3)) – The n three-dimensional coordinates before rotation (used to determine the rectangular box limits).

  • angles (ndarray, shape (2,)) – The two Euler angles for rotation around the first two axes. The last rotation around the normal axis is implicitly set to zero.

  • dir_key (str) – The key indicating the normal Cartesian direction. The rotation around this axis will be the last one and is implicitly set to zero. Must be either "x", "y", or "z".

Returns:

  • B (ndarray, shape (3, 3)) – The vectors spanning the lattice planes before rotation, represented as a matrix with row vectors. The normal vector of the lattice planes is returned with unit length.

  • B_r (ndarray, shape (3, 3)) – The vectors spanning the lattice planes after rotation, represented as a matrix with row vectors. The normal vector of the lattice planes is returned with unit length (and contains only one component in normal direction by definition).

apyt.analysis.crystallography.sdm.optimize_alignment(angles, data, rdf_params, rdf_model_params)[source]

Align lattice planes in normal direction.

For optimal evaluation of the spatial distribution maps, the lattice planes in normal direction should be aligned in an optimal way. This is accomplished with an optimization routine which minimizes the standard deviation of the Gaussian functions as described in rdf(), i.e. the lattice planes are aligned in such a way that the peaks in the radial distribution function (RDF) become the sharpest.

First, the RDF will be generated with the parameters specified in rdf_params (see also rdf_1d()). Then, the RDF will be fitted with the analytical model using the (initial) parameters specified in rdf_model_params (see rdf()).

A simplex algorithm is then used to find the minimal standard deviation by systematically varying the Euler angles around the two lateral axes (see also rotate()).

Parameters:
  • angles (ndarray, shape (2,)) – The initial Euler angles.

  • data (ndarray, shape (n, 3)) – The n three-dimensional coordinates.

  • rdf_params (tuple) – The parameters to generate the RDF as described in rdf_1d() (excluding data).

  • rdf_model_params (tuple) – The parameters used for fitting the RDF as described in rdf() (excluding r). N, n, and w will be fixed and must be exact. The values for sigma and d should be approximate and will be used as initial values for fitting the RDF.

Returns:

angles_opt – The optimal Euler angles.

Return type:

ndarray, shape (2,)

apyt.analysis.crystallography.sdm.rdf(r, N, n, w, sigma, d)[source]

Analytic form of the radial distribution function (RDF).

The RDF is modeled as

\[\mathrm{rdf}(r) = N\, n\, (n - 1)\, w\, f(r, 0, \sigma) + \sum_{i = 1}^{N - 1} (N - i)\, n^2\, w\, f(r, i, d, \sigma),\]

where \(f(r, r_0, \sigma)\) is the Gaussian distribution centered at \(r_0\) with standard deviation \(\sigma\). The fist term represents the contribution from atomic pairs within the same atomic layer and the second term stems from the contributions of atoms in different atomic layers.

Parameters:
  • r (float) – The independent function argument.

  • N (int) – The number of atomic layers in the evaluation volume.

  • n (float) – The number of atoms per layer.

  • w (float) – The width used for the histogram binning (needed for correct normalization).

  • sigma (float) – The standard deviation of the Gaussian distributions.

  • d (float) – The lattice spacing.

Returns:

rdf(r) – The RDF evaluated at r.

Return type:

float

apyt.analysis.crystallography.sdm.rdf_1d(data, min, max, w, dir_key)[source]

Calculate one-dimensional radial distribution function (RDF) for direction specified in dir_key and return histogram data points.

Parameters:
  • data (ndarray, shape (n, 3)) – The n three-dimensional coordinates.

  • min (float) – The minimum distance which should be included in the RDF.

  • max (float) – The maximum distance which should be included in the RDF.

  • w (float) – The width used for binning the histogram data.

  • dir_key (str) – The key indicating the Cartesian direction. Must be either "x", "y", or "z".

Returns:

(x, y) – The touple containing the bin centers x and the respective histogram counts y, both of type ndarray with shape (m,).

Return type:

tuple

apyt.analysis.crystallography.sdm.rdf_lateral(data, min, max, w, dir_key, r_0, δ)[source]

Calculate lateral radial distribution function (RDF) and return histogram data points.

The RDF is calculated for the directions which are orthogonal to the normal direction specified in dir_key, but the lateral RDF is restricted to all pairs for which the distance in normal direction \(r_\mathrm n\) is limited to the range \(r_0 - \frac{\delta}{2} \leq r_\mathrm n \leq r_0 + \frac{\delta}{2}\), i.e. all pairs are chosen which fall in a narrow distance window centered at \(r_0\) in normal direction.

The lateral RDF may reveal peaks for highly precise atomic positions, but this is unlikely for APT data and rather applies to simulation data. The problem with the lateral RDF is that all angular, local in-plane information is averaged out, with the noise dominating the result. The idea behind the SDMs is the exactly the point to keep this local angular information, leading to distinct peaks in the SDMs, even with noisy data.

Parameters:
  • data (ndarray, shape (n, 3)) – The n three-dimensional coordinates.

  • min (float) – The minimum distance which should be included in the RDF.

  • max (float) – The maximum distance which should be included in the RDF.

  • w (float) – The width used for binning the histogram data.

  • dir_key (str) – The key indicating the normal Cartesian direction. Must be either "x", "y", or "z".

  • r_0 (float) – The center of the distance window in normal direction.

  • δ (float) – The width of the distance window in normal direction.

Returns:

(x, y) – The touple containing the bin centers x and the respective histogram counts y, both of type ndarray with shape (m,).

Return type:

tuple

apyt.analysis.crystallography.sdm.rectify_sdm(data, sdm_, dir_key, d, δ, plt_title, use_filter)[source]

Use maxima of spatial distribution map (SDM) for data rectification.

This method asks the user for the (approximate) positions of two linearly independent maxima in the SDM and performs two-dimensional Gaussian fits on these maxima to get accurate values for the centers of the maxima. By comparing these positions to the real crystallographic positions, a system of linear equations is solved to obtain the required transformation matrix which rectifies the SDM (and the atomic positions).

Parameters:
  • data (ndarray, shape (n, 3)) – The n three-dimensional coordinates.

  • sdm (tuple) – The SDM histogram data (H, xedges, yedges), as an element of the list returned by sdm().

  • dir_key (str) – The key indicating the normal Cartesian direction. Must be either "x", "y", or "z".

  • d (float) – The current lattice spacing in normal direction which is required for correct scaling. Should be determined first from fitting the radial distribution function in rdf().

  • δ (list) – The list containing the different widths of the distance window in normal direction, each of type float. Required for the calculation of the rectified SDM(s) (see sdm() for details).

  • plt_title (str) – The title used in generated plots.

  • use_filter (bool) – Whether to apply a filter to the rectified SDMs.

Returns:

  • T (ndarray, shape(3, 3)) – The transformation matrix for data rectification.

  • data_rect (ndarray, shape (n, 3)) – The rectified n three-dimensional coordinates.

  • [sdm_1, sdm_2, …] (list) – The list containing the histogram data (H, xedges, yedges) of all rectified SDMs (see sdm() for details).

apyt.analysis.crystallography.sdm.rotate(data, angles, dir_key)[source]

Rotate three-dimensional data by given Euler angles.

Parameters:
  • data (ndarray, shape (n, 3)) – The n three-dimensional coordinates.

  • angles (ndarray, shape (2,)) – The two Euler angles for rotation around the first two axes. The last rotation around the normal axis is implicitly set to zero.

  • dir_key (str) – The key indicating the normal Cartesian direction. The rotation around this axis will be the last one and is implicitly set to zero. Must be either "x", "y", or "z".

Returns:

data_r – The n three-dimensional coordinates after rotation.

Return type:

ndarray, shape (n, 3)

apyt.analysis.crystallography.sdm.sdm(data, max, n, dir_key, d_0, δ, plt_title, use_filter)[source]

Calculate spatial distribution maps (SDMs) and return the results as a list of two-dimensional histogram data.

The SDMs are calculated for the directions which are orthogonal to the normal direction specified in dir_key, but the (lateral) SDM is restricted to all pairs for which the distance in normal direction \(d\) is limited to the range \(i\, d_0 - \frac{\delta_j}{2} \leq d \leq i\, d_0 + \frac{\delta_j}{2}\), i.e. all pairs are chosen which fall in a narrow distance window centered at \(i\, d_0\) in normal direction. The SDMs are created for \(i \in (0,1,2,3)\), i.e. for pairs in the same layer (\(i = 0\)) to pairs separated by up to three lattice spacings (\(i = 3\)).

Different window widths \(\delta_j\) in normal direction can be provided. For each \(\delta_j\), a common plot is generated and shown for all \(i \in (0,1,2,3)\).

The two-dimensional histogram data of all processed SDMs is returned as a list.

Parameters:
  • data (ndarray, shape (n, 3)) – The n three-dimensional coordinates.

  • max (float) – The maximum distance which should be included in the SDMs.

  • n (int) – The number of bins used for the two-dimensional histogram data.

  • dir_key (str) – The key indicating the normal Cartesian direction. Must be either "x", "y", or "z".

  • d_0 (float) – The lattice spacing in normal direction.

  • δ (list) – The list containing the different widths of the distance window in normal direction, each of type float.

  • plt_title (str) – The plot title.

  • use_filter (bool) – Whether to apply a filter to the SDMs.

Returns:

[sdm_1, sdm_2, …] – The list containing the histogram data (H, xedges, yedges) of all processed SDMs.

Return type:

list