The APyT basic reconstruction module
This module provides methods for performing three-dimensional reconstructions of atom probe tomography (APT) data using the classic and taper geometry models.
It is intended for standard reconstruction scenarios where the tip geometry can be approximated by well-defined analytical models. Both reconstruction approaches rely on geometric assumptions about the specimen shape, field evaporation behavior, and projection parameters.
List of functions
The following functions are provided for tip geometry calculation, evaporation field estimation, and APT data reconstruction:
align_evaporation_field(): Automatically align evaporation field for taper geometry.enable_debug(): Enable or disable debug output.get_evaporation_field(): Calculate evaporation field.get_geometry_classic(): Calculate tip geometry through ‘classic’ scheme.get_geometry_taper(): Calculate taper geometry.reconstruct(): Reconstruct \(xyz\) tip positions.
- apyt.reconstruction.basic.align_evaporation_field(V, U, params_in, E_0, θ_m, num_points=10000)[source]
Automatically align evaporation field for taper geometry.
This function takes \(r_0\), \(θ\), and \(β\) as variation parameters to align the evaporation field to the given target value through least squares minimization.
- Parameters:
V (ndarray, shape (n,)) – The cumulated reconstructed volume for the n events.
U (ndarray, shape (n,)) – The measured voltages of the n events.
params_in (dict) – The parameters to be optimized, represented by the dictionary keys
r_0,θ, andβ. The values are 2-tuples, representing the initial guess (float) and whether the parameter should be varied (bool).E_0 (float) – The target evaporation field for alignment.
θ_m (float) – The (half) aperture angle (in radians).
num_points (int) – The number of points used to sample the evaporation field over the given entire data set. The curve can be considered to change only slowly so that a small subset is sufficient for evaluation. Defaults to
10000sample points.
- Returns:
params_out – The optimized parameters \((r_0, θ, β)\).
- Return type:
tuple
- apyt.reconstruction.basic.enable_debug(is_dbg)[source]
Enable or disable debug output.
- Parameters:
is_dbg (bool) – Whether to enable or disable debug output.
- apyt.reconstruction.basic.get_evaporation_field(U, r, β)[source]
Calculate evaporation field.
The evaporation field is calculated according to the relation \(E = \frac{U}{\beta r}\), where \(\beta\) is the field factor.
- Parameters:
U (ndarray, shape (n,)) – The measured voltages of the n events.
r (ndarray, shape (n,)) – The corresponding radii of curvature.
β (float) – The field factor.
- Returns:
E – The evaporation field of the n events. The precision is inferred from the data type of the U parameter.
- Return type:
ndarray, shape (n,)
- apyt.reconstruction.basic.get_geometry_classic(Ω, r, θ_m)[source]
Calculate tip geometry through ‘classic’ scheme.
This function calculates the tip geometry (current \(z\)-position of sphere center and radius of curvature) through the classical scheme for each event to be reconstructed. The reconstructed volume between two events is calculated exactly with the aid of a local taper geometry according to the two consecutive radii of curvature.
- Parameters:
Ω (ndarray, shape (n,)) – The individual atomic volumes for the n events.
r (ndarray, shape (n,)) – The calculated radii of curvature for the n events.
θ_m (float) – The (half) aperture angle (in radians).
- Returns:
z (ndarray, shape (n,)) – The \(z\)-positions of the sphere centers.
r (ndarray, shape (n,)) – The corresponding radii of curvature.
- apyt.reconstruction.basic.get_geometry_taper(V, r_0, θ, θ_m, use_numba=True)[source]
Calculate taper geometry.
This function solves a cubic equation according to Cardano’s method to calculate the taper geometry (current \(z\)-position of sphere center and radius of curvature) for each event to be reconstructed. Only the physical solution is calculated and returned while the other two complex solutions are dropped. This method is taken in modified form from "Solving Cubic Equation using Cardano’s method with Python" .
- Parameters:
V (ndarray, shape (n,)) – The cumulated reconstructed volume for the n events.
r_0 (float) – The initial radius of the hemispherical cap.
θ (float) – The polar angle of the hemispherical cap (in radians). This angle is related to the taper angle \(\alpha\) through the relation \(\theta = \frac{\pi - \alpha}{2}\).
θ_m (float) – The (half) aperture angle (in radians).
use_numba (bool) – Whether to use the Numba njit’ed function. For low overall computational costs, calling the Numba function may cause too much overhead so that the native Python implementation may be faster, which can be used by setting the use_numba flag to
False. Defaults toTrue, i.e. use the Numba function.
- Returns:
z (ndarray, shape (n,)) – The \(z\)-positions of the sphere centers.
r (ndarray, shape (n,)) – The corresponding radii of curvature.
- apyt.reconstruction.basic.reconstruct(xy_det, z_0, r, L, ξ)[source]
Reconstruct \(xyz\) tip positions.
This function reconstructs the \(xyz\) tip positions using the conventional point projection algorithm from Bas et al.
- Parameters:
xy_det (ndarray, shape (n, 2)) – The \(xy\) detector positions for the n events.
z_0 (ndarray, shape (n,)) – The \(z\)-positions of the sphere centers.
r (ndarray, shape (n,)) – The corresponding radii of curvature.
L (float) – The distance between the tip and the detector. The unit must match the unit of the \(xy\) detector positions.
ξ (float) – The image compression factor which is defined as \(\xi = \frac{\theta_{\textnormal{crys}}}{\tan\theta_{\textnormal{obs}}}\).
- Returns:
xyz_tip – The \(xyz\) tip positions of the n reconstructed events. The precision is inferred from the data type of the xy_det parameter.
- Return type:
ndarray, shape (n, 3)