abel package

abel.transform module

class abel.transform.Transform(IM, direction=u'inverse', method=u'three_point', center=u'none', symmetry_axis=None, use_quadrants=(True, True, True, True), symmetrize_method=u'average', angular_integration=False, transform_options={}, center_options={}, angular_integration_options={}, recast_as_float64=True, verbose=False)[source]

Bases: object

Abel transform image class.

This class provides whole image forward and inverse Abel transformations, together with preprocessing (centering, symmetrizing) and post processing (integration) functions.

The following class attributes are available, depending on the calculation.

Returns:
  • transform (numpy 2D array) – the 2D forward/reverse Abel transform.
  • angular_integration (tuple) – (radial-grid, radial-intensity) radial coordinates, and the radial intensity (speed) distribution, evaluated using abel.tools.vmi.angular_integration().
  • residual (numpy 2D array) – residual image (not currently implemented).
  • IM (numpy 2D array) – the input image, re-centered (optional) with an odd-size width.
  • method (str) – transform method, as specified by the input option.
  • direction (str) – transform direction, as specified by the input option.
  • Beta (numpy 2D array) – with linbasex transform_options=dict(return_Beta=True)() Beta array coefficients of Newton sphere spherical harmonics
    Beta[0] - the radial intensity variation

    Beta[1] - the anisotropy parameter variation

    …Beta[n] - higher order terms up to legedre_orders = [0, …, n]

  • radial (numpy 1d array) – with linbasex transform_options=dict(return_Beta=True)() radial-grid for Beta array
  • projection – with linbasex transform_options=dict(return_Beta=True)() radial projection profiles at angles proj_angles
__init__(IM, direction=u'inverse', method=u'three_point', center=u'none', symmetry_axis=None, use_quadrants=(True, True, True, True), symmetrize_method=u'average', angular_integration=False, transform_options={}, center_options={}, angular_integration_options={}, recast_as_float64=True, verbose=False)[source]

The one stop transform function.

Parameters:
  • IM (a NxM numpy array) – This is the image to be transformed
  • direction (str) –

    The type of Abel transform to be performed.

    forward
    A ‘forward’ Abel transform takes a (2D) slice of a 3D image and returns the 2D projection.
    inverse
    An ‘inverse’ Abel transform takes a 2D projection and reconstructs a 2D slice of the 3D image.

    The default is inverse.

  • method (str) –

    specifies which numerical approximation to the Abel transform should be employed (see below). The options are

    hansenlaw
    the recursive algorithm described by Hansen and Law.
    basex
    the Gaussian “basis set expansion” method of Dribinski et al.
    direct
    a naive implementation of the analytical formula by Roman Yurchuk.
    two_point
    the two-point transform of Dasch (1992).
    three_point
    the three-point transform of Dasch (1992).
    onion_bordas
    the algorithm of Bordas and co-workers (1996), re-implemented by Rallis, Wells and co-workers (2014).
    onion_peeling
    the onion peeling deconvolution as described by Dasch (1992).
    linbasex
    the 1d-projections of VM-images in terms of 1d spherical functions by Gerber et al. (2013).
  • center (tuple or str) –

    If a tuple (float, float) is provided, this specifies the image center in (y,x) (row, column) format. A value None can be supplied if no centering is desired in one dimension, for example ‘center=(None, 250)’. If a string is provided, an automatic centering algorithm is used

    image_center
    center is assumed to be the center of the image.
    convolution
    center the image by convolution of two projections along each axis.
    slice
    the center is found my comparing slices in the horizontal and vertical directions
    com
    the center is calculated as the center of mass
    gaussian
    the center is found using a fit to a Gaussian function. This only makes sense if your data looks like a Gaussian.
    none
    (Default) No centering is performed. An image with an odd number of columns must be provided.
  • symmetry_axis (None, int or tuple) – Symmetrize the image about the numpy axis 0 (vertical), 1 (horizontal), (0,1) (both axes)
  • use_quadrants (tuple of 4 booleans) – select quadrants to be used in the analysis: (Q0,Q1,Q2,Q3). Quadrants are numbered counter-clockwide from upper right. See note below for description of quadrants. Default is (True, True, True, True), which uses all quadrants.
  • symmetrize_method (str) –

    Method used for symmetrizing the image.

    average
    average the quadrants, in accordance with the symmetry_axis
    fourier
    axial symmetry implies that the Fourier components of the 2-D projection should be real. Removing the imaginary components in reciprocal space leaves a symmetric projection. ref: Overstreet, K., et al. “Multiple scattering and the density distribution of a Cs MOT.” Optics express 13.24 (2005): 9672-9682. http://dx.doi.org/10.1364/OPEX.13.009672
  • angular_integration (boolean) – integrate the image over angle to give the radial (speed) intensity distribution
  • transform_options (tuple) – Additional arguments passed to the individual transform functions. See the documentation for the individual transform method for options.
  • center_options (tuple) – Additional arguments to be passed to the centering function.
  • angular_integration_options (tuple (or dict)) – Additional arguments passed to the angular_integration transform functions. See the documentation for angular_integration for options.
  • recast_as_float64 (boolean) – True/False that determines if the input image should be recast to float64. Many images are imported in other formats (such as uint8 or uint16) and this does not always play well with the transorm algorithms. This should probably always be set to True. (Default is True.)
  • verbose (boolean) – True/False to determine if non-critical output should be printed.

Note

Quadrant combining

The quadrants can be combined (averaged) using the use_quadrants keyword in order to provide better data quality.

The quadrants are numbered starting from Q0 in the upper right and proceeding counter-clockwise:

+--------+--------+
| Q1   * | *   Q0 |
|   *    |    *   |
|  *     |     *  |                               AQ1 | AQ0
+--------o--------+ --(inverse Abel transform)--> ----o----
|  *     |     *  |                               AQ2 | AQ3
|   *    |    *   |
| Q2  *  | *   Q3 |          AQi == inverse Abel transform
+--------+--------+                 of quadrant Qi

Three cases are possible:

  1. symmetry_axis = 0 (vertical):

    Combine:  Q01 = Q0 + Q1, Q23 = Q2 + Q3
    inverse image   AQ01 | AQ01
                    -----o----- (left and right sides equivalent)
                    AQ23 | AQ23
    
  2. symmetry_axis = 1 (horizontal):

    Combine: Q12 = Q1 + Q2, Q03 = Q0 + Q3
    inverse image   AQ12 | AQ03
                    -----o----- (top and bottom equivalent)
                    AQ12 | AQ03
    
  3. symmetry_axis = (0, 1) (both):

    Combine: Q = Q0 + Q1 + Q2 + Q3
    inverse image   AQ | AQ
                    ---o---  (all quadrants equivalent)
                    AQ | AQ
    

Notes

As mentioned above, PyAbel offers several different approximations to the the exact abel transform. All the the methods should produce similar results, but depending on the level and type of noise found in the image, certain methods may perform better than others. Please see the “Transform Methods” section of the documentation for complete information.

hansenlaw

This “recursive algorithm” produces reliable results and is quite fast (~0.1 sec for a 1001x1001 image). It makes no assumptions about the data (apart from cylindrical symmetry). It tends to require that the data is finely sampled for good convergence.

E. W. Hansen and P.-L. Law “Recursive methods for computing the Abel transform and its inverse” J. Opt. Soc. A*2, 510-520 (1985) http://dx.doi.org/10.1364/JOSAA.2.000510

basex *

The “basis set exapansion” algorithm describes the data in terms of gaussian functions, which themselves can be abel transformed analytically. Because the gaussian functions are approximately the size of each pixel, this method also does not make any assumption about the shape of the data. This method is one of the de-facto standards in photoelectron/photoion imaging.

Dribinski et al, 2002 (Rev. Sci. Instrum. 73, 2634) http://dx.doi.org/10.1063/1.1482156
direct
This method attempts a direct integration of the Abel transform integral. It makes no assumptions about the data (apart from cylindrical symmetry), but it typically requires fine sampling to converge. Such methods are typically inefficient, but thanks to this Cython implementation (by Roman Yurchuk), this ‘direct’ method is competitive with the other methods.
linbasex *

VM-images are composed of projected Newton spheres with a common centre. The 2D images are usually evaluated by a decomposition into base vectors each representing the 2D projection of a set of particles starting from a centre with a specific velocity distribution. linbasex evaluate 1D projections of VM-images in terms of 1D projections of spherical functions, instead.

..Rev. Sci. Instrum. 84, 033101 (2013): <http://scitation.aip.org/content/aip/journal/rsi/84/3/10.1063/1.4793404>

onion_bordas
The onion peeling method, also known as “back projection”, originates from Bordas et al. Rev. Sci. Instrum. 67, 2257 (1996).
The algorithm was subsequently coded in MatLab by Rallis, Wells and co-workers, Rev. Sci. Instrum. 85, 113105 (2014).
which was used as the basis of this Python port. See issue #56.
onion_peeling *
This is one of the most compact and fast algorithms, with the inverse Abel transfrom achieved in one Python code-line, PR #155. See also three_point is the onion peeling algorithm as described by Dasch (1992), reference below.
two_point *
Another Dasch method. Simple, and fast, but not as accurate as the other methods.
three_point *

The “Three Point” Abel transform method exploits the observation that the value of the Abel inverted data at any radial position r is primarily determined from changes in the projection data in the neighborhood of r. This method is also very efficient once it has generated the basis sets.

Dasch, 1992 (Applied Optics, Vol 31, No 8, March 1992, Pg 1146-1152).

*
The methods marked with a * indicate methods that generate basis sets. The first time they are run for a new image size, it takes seconds to minutes to generate the basis set. However, this basis set is saved to disk can can be reloaded, meaning that future transforms are performed much more quickly.
__weakref__

list of weak references to the object (if defined)

abel.basex module

abel.basex.basex_transform(data, nbf=u'auto', basis_dir=u'./', dr=1.0, verbose=True, direction=u'inverse')[source]

This function performs the BASEX (BAsis Set EXpansion) Abel Transform. It works on a “right side” image. I.e., it works on just half of a cylindrically symmetric object and data[0,0] should correspond to a central pixel. To perform a BASEX transorm on a whole image, use

abel.Transform(image, method='basex', direction='inverse').transform

This BASEX implementation only works with images that have an odd-integer width.

Note: only the direction=”inverse” transform is currently implemented.

Parameters:
  • data (a NxM numpy array) – The image to be inverse transformed. The width (M) must be odd and data[:,0] should correspond to the central column of the image.
  • nbf (str or list) – number of basis functions. If nbf='auto', it is set to (n//2 + 1). This is what you should always use, since this BASEX implementation does not work reliably in other situations! In the future, you could use list, in format [nbf_vert, nbf_horz]
  • basis_dir (str) – path to the directory for saving / loading the basis set coefficients. If None, the basis set will not be saved to disk.
  • dr (float) – size of one pixel in the radial direction. This only affects the absolute scaling of the transformed image.
  • verbose (boolean) – Determines if statements should be printed.
  • direction (str) – The type of Abel transform to be performed. Currently only accepts value 'inverse'.
Returns:

recon – the transformed (half) image

Return type:

NxM numpy array

abel.basex.basex_core_transform(rawdata, M_vert, M_horz, Mc_vert, Mc_horz, vert_left, horz_right, dr=1.0)[source]

This is the internal function that does the actual BASEX transform. It requires that the matrices of basis set coefficients be passed.

Parameters:
  • rawdata (NxM numpy array) – the raw image. This is the full image, both left and right sides.
  • M_vert_etc. (Numpy arrays) – 2D arrays given by the basis set calculation function
  • dr (float) – pixel size. This only affects the absolute scaling of the output.
Returns:

IM – The abel-transformed image, a slice of the 3D distribution

Return type:

NxM numpy array

abel.basex.get_bs_basex_cached(n_vert, n_horz, nbf=u'auto', basis_dir=u'.', verbose=False)[source]

Internal function.

Gets BASEX basis sets, using the disk as a cache (i.e. load from disk if they exist, if not calculate them and save a copy on disk). To prevent saving the basis sets to disk, set basis_dir=None.

Parameters:
  • n_horz (n_vert,) – Abel inverse transform will be performed on a n_vert x n_horz area of the image
  • nbf (int or list) – number of basis functions. If nbf='auto', n_horz is set to (n//2 + 1).
  • basis_dir (str) – path to the directory for saving / loading the basis set coefficients. If None, the basis sets will not be saved to disk.
Returns:

M_vert, M_horz, Mc_vert, Mc_horz, vert_left, horz_right – the matrices that compose the basis set.

Return type:

numpy arrays

abel.linbasex module

abel.linbasex.linbasex_transform(IM, proj_angles=[0, 1.5707963267948966], legendre_orders=[0, 2], radial_step=1, smoothing=0.5, rcond=0.0005, threshold=0.2, basis_dir=None, return_Beta=False, clip=0, norm_range=(0, -1), direction=u'inverse', verbose=False, **kwargs)[source]

wrapper function for linebasex to process supplied quadrant-image as a full-image.

PyAbel transform functions operate on the right side of an image. Here we follow the basex technique of duplicating the right side to the left re-forming the whole image.

Inverse Abel transform using 1d projections of images.

Gerber, Thomas, Yuzhu Liu, Gregor Knopp, Patrick Hemberger, Andras Bodi, Peter Radi, and Yaroslav Sych. Charged Particle Velocity Map Image Reconstruction with One-Dimensional Projections of Spherical Functions. Review of Scientific Instruments 84, no. 3 (March 1, 2013): 033101–033101 – 10. <http://dx.doi.org/10.1063/1.4793404>`_

``linbasex``models the image using a sum of Legendre polynomials at each radial pixel, As such, it should only be applied to situations that can be adequately represented by Legendre polynomials, i.e., images that feature spherical-like structures. The reconstructed 3D object is obtained by adding all the contributions, from which slices are derived.
Parameters:
  • IM (numpy 2D array) – image data must be square shape of odd size
  • proj_angles (list) – projection angles, in radians (default \([0, \pi/2]\)) e.g. \([0, \pi/2]\) or \([0, 0.955, \pi/2]\) or \([0, \pi/4, \pi/2, 3\pi/4]\)
  • legendre_orders (list) –
    orders of Legendre polynomials to be used as the expansion
    even polynomials [0, 2, …] gerade odd polynomials [1, 3, …] ungerade all orders [0, 1, 2, …].

    In a single photon experiment there are only anisotropies up to second order. The interaction of 4 photons (four wave mixing) yields anisotropies up to order 8.

  • radial_step (int) – number of pixels per Newton sphere (default 1)
  • smoothing (float) – convolve Beta array with a Gaussian function of 1/e 1/2 width smoothing.
  • rcond (float) – (default 0.0005) scipy.linalg.lstsq fit conditioning value. set rcond to zero to switch conditioning off. Note: In the presence of noise the equation system may be ill posed. Increasing rcond smoothes the result, lowering it beyond a minimum renders the solution unstable. Tweak rcond to get a “reasonable” solution with acceptable resolution.
  • clip (int) – clip first vectors (smallest Newton spheres) to avoid singularities (default 0)
  • norm_range (tuple) – (low, high) normalization of Newton spheres, maximum in range Beta[0, low:high]. Note: Beta[0, i] the total number of counts integrated over sphere i, becomes 1.
  • threshold (float) – threshold for normalization of higher order Newton spheres (default 0.2) Set all Beta[j], j>=1 to zero if the associated Beta[0] is smaller than threshold.
  • return_Beta (bool) –

    return the Beta array of Newton spheres, as the tuple: radial-grid, Beta for the case legendre_orders=[0, 2]

    Beta[0] vs radius -> speed distribution

    Beta[2] vs radius -> anisotropy of each Newton sphere

    see ‘Returns’.

  • direction (str) – “inverse” - only option for this method. Abel transform direction.
  • verbose (bool) – print information about processing (normally used for debugging)
Returns:

  • inv_IM (numpy 2D array) – inverse Abel transformed image

  • radial-grid, Beta, projections (tuple) – (if return_Beta=True)

    contributions of each spherical harmonic \(Y_{i0}\) to the 3D distribution contain all the information one can get from an experiment. For the case legendre_orders=[0, 2]:

    Beta[0] vs radius -> speed distribution

    Beta[1] vs radius -> anisotropy of each Newton sphere.

    projections : are the radial projection profiles at angles proj_angles

abel.linbasex.linbasex_transform_full(IM, proj_angles=[0, 1.5707963267948966], legendre_orders=[0, 2], radial_step=1, smoothing=0.5, rcond=0.0005, threshold=0.2, clip=0, basis_dir=None, return_Beta=False, norm_range=(0, -1), direction=u'inverse', verbose=False, **kwargs)[source]

Inverse Abel transform using 1d projections of images.

Gerber, Thomas, Yuzhu Liu, Gregor Knopp, Patrick Hemberger, Andras Bodi, Peter Radi, and Yaroslav Sych. Charged Particle Velocity Map Image Reconstruction with One-Dimensional Projections of Spherical Functions. Review of Scientific Instruments 84, no. 3 (March 1, 2013): 033101–033101 – 10. <http://dx.doi.org/10.1063/1.4793404>`_

``linbasex``models the image using a sum of Legendre polynomials at each radial pixel, As such, it should only be applied to situations that can be adequately represented by Legendre polynomials, i.e., images that feature spherical-like structures. The reconstructed 3D object is obtained by adding all the contributions, from which slices are derived.
Parameters:
  • IM (numpy 2D array) – image data must be square shape of odd size
  • proj_angles (list) – projection angles, in radians (default \([0, \pi/2]\)) e.g. \([0, \pi/2]\) or \([0, 0.955, \pi/2]\) or \([0, \pi/4, \pi/2, 3\pi/4]\)
  • legendre_orders (list) –
    orders of Legendre polynomials to be used as the expansion
    even polynomials [0, 2, …] gerade odd polynomials [1, 3, …] ungerade all orders [0, 1, 2, …].

    In a single photon experiment there are only anisotropies up to second order. The interaction of 4 photons (four wave mixing) yields anisotropies up to order 8.

  • radial_step (int) – number of pixels per Newton sphere (default 1)
  • smoothing (float) – convolve Beta array with a Gaussian function of 1/e 1/2 width smoothing.
  • rcond (float) – (default 0.0005) scipy.linalg.lstsq fit conditioning value. set rcond to zero to switch conditioning off. Note: In the presence of noise the equation system may be ill posed. Increasing rcond smoothes the result, lowering it beyond a minimum renders the solution unstable. Tweak rcond to get a “reasonable” solution with acceptable resolution.
  • clip (int) – clip first vectors (smallest Newton spheres) to avoid singularities (default 0)
  • norm_range (tuple) – (low, high) normalization of Newton spheres, maximum in range Beta[0, low:high]. Note: Beta[0, i] the total number of counts integrated over sphere i, becomes 1.
  • threshold (float) – threshold for normalization of higher order Newton spheres (default 0.2) Set all Beta[j], j>=1 to zero if the associated Beta[0] is smaller than threshold.
  • return_Beta (bool) –

    return the Beta array of Newton spheres, as the tuple: radial-grid, Beta for the case legendre_orders=[0, 2]

    Beta[0] vs radius -> speed distribution

    Beta[2] vs radius -> anisotropy of each Newton sphere

    see ‘Returns’.

  • direction (str) – “inverse” - only option for this method. Abel transform direction.
  • verbose (bool) – print information about processing (normally used for debugging)
Returns:

  • inv_IM (numpy 2D array) – inverse Abel transformed image

  • radial-grid, Beta, projections (tuple) – (if return_Beta=True)

    contributions of each spherical harmonic \(Y_{i0}\) to the 3D distribution contain all the information one can get from an experiment. For the case legendre_orders=[0, 2]:

    Beta[0] vs radius -> speed distribution

    Beta[1] vs radius -> anisotropy of each Newton sphere.

    projections : are the radial projection profiles at angles proj_angles

abel.linbasex.int_beta(Beta, radial_step=1, threshold=0.1, regions=None)[source]

Integrate beta over a range of Newton spheres.

Parameters:
  • Beta (numpy array) – Newton spheres
  • radial_step (int) – number of pixels per Newton sphere (default 1)
  • threshold (float) – threshold for normalisation of higher orders, 0.0 … 1.0.
  • regions (list of tuple radial ranges) – [(min0, max0), (min1, max1), …]
Returns:

Beta_in – integrated normalized Beta array [Newton sphere, region]

Return type:

numpy array

abel.hansenlaw module

abel.hansenlaw.hansenlaw_transform(image, dr=1, direction=u'inverse', hold_order=0, sub_pixel_shift=-0.35, **kwargs)[source]

Forward/Inverse Abel transformation using the algorithm of:

E. W. Hansen “Fast Hankel Transform” IEEE Trans. Acoust. Speech Signal Proc. 33, 666 (1985)

and

E. W. Hansen and P.-L. Law “Recursive methods for computing the Abel transform and its inverse” J. Opt. Soc. Am. A 2, 510-520 (1985)

This function performs the Hansen-Law transform on only one “right-side” image:

Abeltrans = abel.hansenlaw.hansenlaw_transform(image, direction='inverse')

Note

Image should be a right-side image, like this:

.         +--------      +--------+
.         |      *       | *      |
.         |   *          |    *   |  <---------- im
.         |  *           |     *  |
.         +--------      o--------+
.         |  *           |     *  |
.         |   *          |    *   |
.         |     *        | *      |
.         +--------      +--------+

In accordance with all PyAbel methods the image center o is defined to be mid-pixel i.e. an odd number of columns, for the full image.

For the full image transform, use the :class:abel.Transform.

Inverse Abel transform:

iAbel = abel.Transform(image, method='hansenlaw').transform

Forward Abel transform:

fAbel = abel.Transform(image, direction='forward', method='hansenlaw').transform
Parameters:
  • image (1D or 2D numpy array) – Right-side half-image (or quadrant). See figure below.
  • dr (float) – Sampling size, used for Jacobian scaling. Default: 1 (appliable for pixel images).
  • direction (string 'forward' or 'inverse') – forward or inverse Abel transform. Default: ‘inverse’.
  • hold_order (int 0 or 1) – The order of the hold approximation, used to evaluate the state equation integral. 0 assumes a constant intensity across a pixel (between grid points) for the driving function (the image gradient for the inverse transform, or the original image, for the forward transform). 1 assumes a linear intensity variation between grid points, which may yield a more accurate transform for some functions (see PR 211). Default: 0.
  • sub_pixel_shift (float) – For the zero-order hold approximation hold_order=0, a sub-pixel left-shift of the driving function (image-forward or gradient-inverse) improves the transform alignment with the other PyAbel methods, and Abel transform-pair functions. See the discussion in final summary of PR #211 Default: -0.35.
Returns:

aim – forward/inverse Abel transform half-image

Return type:

1D or 2D numpy array

abel.dasch module

abel.dasch.two_point_transform(IM, basis_dir=u'.', dr=1, direction=u'inverse')[source]
two-point deconvolution
C. J. Dasch Applied Optics 31, 1146 (1992). http://dx.doi.org/10.1364/AO.31.001146
Parameters:
  • IM (1D or 2D numpy array) – right-side half-image (or quadrant)
  • basis_dir (str) – path to the directory for saving / loading the “two-point” operator matrix. If None, the operator matrix will not be saved to disk.
  • dr (float) – sampling size (=1 for pixel images), used for Jacobian scaling. The resulting inverse transform is simply scaled by 1/dr.
  • direction (str) – only the direction=”inverse” transform is currently implemented
Returns:

inv_IM – the “two-point” inverse Abel transformed half-image

Return type:

1D or 2D numpy array

abel.dasch.three_point_transform(IM, basis_dir=u'.', dr=1, direction=u'inverse')[source]
three-point deconvolution
C. J. Dasch Applied Optics 31, 1146 (1992). http://dx.doi.org/10.1364/AO.31.001146
Parameters:
  • IM (1D or 2D numpy array) – right-side half-image (or quadrant)
  • basis_dir (str) – path to the directory for saving / loading the “three-point” operator matrix. If None, the operator matrix will not be saved to disk.
  • dr (float) – sampling size (=1 for pixel images), used for Jacobian scaling. The resulting inverse transform is simply scaled by 1/dr.
  • direction (str) – only the direction=”inverse” transform is currently implemented
Returns:

inv_IM – the “three-point” inverse Abel transformed half-image

Return type:

1D or 2D numpy array

abel.dasch.onion_peeling_transform(IM, basis_dir=u'.', dr=1, direction=u'inverse')[source]
onion-peeling deconvolution
C. J. Dasch Applied Optics 31, 1146 (1992). http://dx.doi.org/10.1364/AO.31.001146
Parameters:
  • IM (1D or 2D numpy array) – right-side half-image (or quadrant)
  • basis_dir (str) – path to the directory for saving / loading the “onion-peeling” operator matrix. If None, the operator matrix will not be saved to disk.
  • dr (float) – sampling size (=1 for pixel images), used for Jacobian scaling. The resulting inverse transform is simply scaled by 1/dr.
  • direction (str) – only the direction=”inverse” transform is currently implemented
Returns:

inv_IM – the “onion-peeling” inverse Abel transformed half-image

Return type:

1D or 2D numpy array

abel.dasch.dasch_transform(IM, D)[source]

Inverse Abel transform using a given D-operator basis matrix.

Parameters:
  • IM (2D numpy array) – image data
  • D (2D numpy array) – D-operator basis shape (cols, cols)
Returns:

inv_IM – inverse Abel transform according to basis operator D

Return type:

2D numpy array

abel.onion_bordas module

abel.onion_bordas.onion_bordas_transform(IM, dr=1, direction=u'inverse', shift_grid=False, **kwargs)[source]

Onion peeling (or back projection) inverse Abel transform.

This algorithm was adapted by Dan Hickstein from the original Matlab implementation, created by Chris Rallis and Eric Wells of Augustana University, and described in this paper:

http://scitation.aip.org/content/aip/journal/rsi/85/11/10.1063/1.4899267

The algorithm actually originates from this 1996 RSI paper by Bordas et al:

http://scitation.aip.org/content/aip/journal/rsi/67/6/10.1063/1.1147044

This function operates on the “right side” of an image. i.e. it works on just half of a cylindrically symmetric image. Unlike the other transforms, the left edge should be the image center, not mid-first pixel. This corresponds to an even-width full image. If not, set shift_grid=True.

To perform a onion-peeling transorm on a whole image, use

abel.Transform(image, method='onion_bordas').transform
Parameters:
  • IM (1D or 2D numpy array) – right-side half-image (or quadrant)
  • dr (float) – sampling size (=1 for pixel images), used for Jacobian scaling. The resulting inverse transform is simply scaled by 1/dr.
  • direction (str) – only the direction=”inverse” transform is currently implemented
  • shift_grid (boolean) – place width-center on grid (bottom left pixel) by shifting image center (-1/2, -1/2) pixel
Returns:

AIM – the inverse Abel transformed half-image

Return type:

1D or 2D numpy array

abel.direct module

abel.direct.direct_transform(fr, dr=None, r=None, direction=u'inverse', derivative=<function gradient>, int_func=<function trapz>, correction=True, backend=u'C', **kwargs)[source]

This algorithm performs a direct computation of the Abel transform integrals. When correction=False, the pixel at the lower bound of the integral (where y=r) is skipped, which causes a systematic error in the Abel transform. However, if correction=True is used, then an analytical transform transform is applied to this pixel, which makes the approximation that the function is linear across this pixel. With correction=True, the Direct method produces reasonable results.

The Direct method is implemented in both Python and, if Cython is available during PyAbel’s installation, a compiled C version, which is much faster. The implementation can be selected using the backend argument.

By default, integration at all other pixels is performed using the Trapezoidal rule.

Parameters:
  • fr (1d or 2d numpy array) – input array to which direct/inverse Abel transform will be applied. For a 2d array, the first dimension is assumed to be the z axis and the second the r axis.
  • dr (float) – spatial mesh resolution (optional, default to 1.0)
  • r (1D ndarray) – the spatial mesh (optional). Unusually, direct_transform should, in principle, be able to handle non-uniform data. However, this has not been regorously tested.
  • direction (string) – Determines if a forward or inverse Abel transform will be applied. can be ‘forward’ or ‘inverse’.
  • derivative (callable) – a function that can return the derivative of the fr array with respect to r. (only used in the inverse Abel transform).
  • int_func (function) – This function is used to complete the integration. It should resemble np.trapz, in that it must be callable using axis=, x=, and dx= keyword arguments.
  • correction (boolean) – If False the pixel where the weighting function has a singular value (where r==y) is simply skipped, causing a systematic under-estimation of the Abel transform. If True, integration near the singular value is performed analytically, by assuming that the data is linear across that pixel. The accuracy of this approximation will depend on how the data is sampled.
  • backend (string) – There are currently two implementations of the Direct transform, one in pure Python and one in Cython. The backend paremeter selects which method is used. The Cython code is converted to C and compiled, so this is faster. Can be ‘C’ or ‘python’ (case insensitive). ‘C’ is the default, but ‘python’ will be used if the C-library is not available.
Returns:

out – with either the direct or the inverse abel transform.

Return type:

1d or 2d numpy array of the same shape as fr

abel.direct.is_uniform_sampling(r)[source]

Returns True if the array is uniformly spaced to within 1e-13. Otherwise False.

Image processing tools

abel.tools.analytical module

class abel.tools.analytical.BaseAnalytical(n, r_max, symmetric=True, **args)[source]

Bases: object

class abel.tools.analytical.StepAnalytical(n, r_max, r1, r2, A0=1.0, ratio_valid_step=1.0, symmetric=True)[source]

Bases: abel.tools.analytical.BaseAnalytical

abel_step_analytical(r, A0, r0, r1)[source]

Directed Abel transform of a step function located between r0 and r1, with a height A0

A0 +                  +-------------+
   |                  |             |
   |                  |             |
 0 | -----------------+             +-------------
   +------------------+-------------+------------>
   0                  r0            r1           r axis
Parameters:
  • r1 (r0,) – vecor of positions along the r axis. Must start with 0.
  • r1 – positions of the step along the r axis
  • A0 (float or 1D array:) – height of the step. If 1D array the height can be variable along the Z axis
Returns:

Return type:

1D array if A0 is a float, a 2D array otherwise

sym_abel_step_1d(r, A0, r0, r1)[source]

Produces a symmetrical analytical transform of a 1d step

class abel.tools.analytical.GaussianAnalytical(n, r_max, sigma=1.0, A0=1.0, ratio_valid_sigma=2.0, symmetric=True)[source]

Bases: abel.tools.analytical.BaseAnalytical

class abel.tools.analytical.TransformPair(n, profile=5)[source]

Bases: abel.tools.analytical.BaseAnalytical

Abel transform pair analytical functions.

profiles1-7: Table 1 of Chan and Hieftje Spectrochimica Acta B 61, 31-41 (2006)

profile8: curve B Hansen and Law J. Opt. Soc. Am. A 2, 510-520 (1985)

Returns:
  • r (numpy array) – vector of positions along the r axis: linspace(0, 1, n)
  • dr (float) – radial interval
  • func (numpy array) – values of the original function (same shape as r)
  • abel (numpy array) – values of the Abel transform (same shape as func)
  • label (str) – name of the curve
  • mask_valid (boolean array) – set all True. Used for unit tests
class abel.tools.analytical.SampleImage(n=361, name='dribinski', sigma=2, temperature=200)[source]

Bases: abel.tools.analytical.BaseAnalytical

abel.tools.basis module

abel.tools.basis.get_bs_cached(method, cols, basis_dir=u'.', basis_options={}, verbose=False)[source]

load basis set from disk, generate and store if not available.

Checks whether file {method}_basis_{cols}_{cols}*.npy is present in basis_dir (*) special case for linbasex

_{legendre_orders}_{proj_angles}_{radial_step}_{clip}
where {legendre_orders} = str of the list elements, typically ‘02’
(proj_angles} = str of the list elements, typically ‘04590135’ {radial_step} = pixel grid size, usually 1 {clip} = clipping size, usually 0

Either, read basis array or generate basis, saving it to the file.

Parameters:
  • method (str) – Abel transform method, currently linbasex, onion_peeling, three_point, and two_point
  • cols (int) – width of image
  • basis_dir (str) – path to the directory for saving / loading the basis
  • verbose (boolean) – print information for debugging
Returns:

  • D (numpy 2D array of shape (cols, cols)) – basis operator array
  • file.npy (file) – saves basis to file name {method}_basis_{cols}_{cols}*.npy * == __{legendre_orders}_{proj_angles}_{radial_step}_{clip} for linbasex method

abel.tools.center module

abel.tools.center.find_center(IM, center=u'image_center', square=False, verbose=False, **kwargs)[source]
Find the coordinates of image center, using the method
specified by the center parameter.
Parameters:
  • IM (2D np.array) – image data
  • center (str) –

    this determines how the center should be found. The options are:

    image_center
    the center of the image is used as the center. The trivial result.
    com
    the center is found as the center of mass.
    convolution
    center by convolution of two projections along each axis.
    gaussian
    the center is extracted by a fit to a Gaussian function. This is probably only appropriate if the data resembles a gaussian.
    slice
    the image is broken into slices, and these slices compared for symmetry.
  • square (bool) – if ‘True’ returned image will have a square shape
Returns:

out – coordinate of the center of the image in the (y,x) format (row, column)

Return type:

(float, float)

abel.tools.center.center_image(IM, center=u'com', odd_size=True, square=False, axes=(0, 1), crop=u'maintain_size', verbose=False, **kwargs)[source]
Center image with the custom value or by several methods provided in
find_center function
Parameters:
  • IM (2D np.array) – The image data.
  • center (tuple or str) –

    center can either be (float, float), the coordinate of the center of the image in the (y,x) format (row, column)

    Or, it can be a string, to specify an automatic centering method. The options are:

    image_center
    the center of the image is used as the center. The trivial result.
    com
    the center is found as the center of mass.
    convolution
    center by convolution of two projections along each axis.
    gaussian
    the center is extracted by a fit to a Gaussian function. This is probably only appropriate if the data resembles a gaussian.
    slice
    the image is broken into slices, and these slices compared for symmetry.
  • odd_size (boolean) – if True, an image will be returned containing an odd number of columns. Most of the transform methods require this, so it’s best to set this to True if the image will subsequently be Abel transformed.
  • square (bool) – if ‘True’ returned image will have a square shape
  • crop (str) –

    This determines how the image should be cropped. The options are:

    maintain_size
    return image of the same size. Some of the original image may be lost and some regions may be filled with zeros.
    valid_region
    return the largest image that can be created without padding. All of the returned image will correspond to the original image. However, portions of the original image will be lost. If you can tolerate clipping the edges of the image, this is probably the method to choose.
    maintain_data
    the image will be padded with zeros such that none of the original image will be cropped.
  • axes (int or tuple) – center image with respect to axis 0 (vertical), 1 (horizontal), or both axes (0, 1) (default).
Returns:

out – Centered image

Return type:

2D np.array

abel.tools.center.set_center(data, center, crop=u'maintain_size', axes=(0, 1), verbose=False)[source]

Move image center to mid-point of image

Parameters:
  • data (2D np.array) – The image data
  • center (tuple) – image pixel coordinate center (row, col)
  • crop (str) –

    This determines how the image should be cropped. The options are:

    maintain_size
    return image of the same size. Some of the original image may be lost and some regions may be filled with zeros.
    valid_region
    return the largest image that can be created without padding. All of the returned image will correspond to the original image. However, portions of the original image will be lost. If you can tolerate clipping the edges of the image, this is probably the method to choose.
    maintain_data
    the image will be padded with zeros such that none of the original image will be cropped.
  • axes (int or tuple) – center image with respect to axis 0 (vertical), 1 (horizontal), or both axes (0, 1) (default).
  • verbose (boolean) – True: print diagnostics
abel.tools.center.find_center_by_center_of_mass(IM, verbose=False, round_output=False, **kwargs)[source]

Find image center by calculating its center of mass

abel.tools.center.find_center_by_convolution(IM, **kwargs)[source]
Center the image by convolution of two projections along each axis.
Code from the linbasex juptyer notebook
IM: numpy 2D array
image data
Returns:center – (row-center, col-center)
Return type:tuple
abel.tools.center.find_center_by_center_of_image(data, verbose=False, **kwargs)[source]

Find image center simply from its dimensions.

abel.tools.center.find_center_by_gaussian_fit(IM, verbose=False, round_output=False, **kwargs)[source]

Find image center by fitting the summation along x and y axis of the data to two 1D Gaussian function.

abel.tools.center.axis_slices(IM, radial_range=(0, -1), slice_width=10)[source]

returns vertical and horizontal slice profiles, summed across slice_width.

Parameters:
  • IM (2D np.array) – image data
  • radial_range (tuple floats) – (rmin, rmax) range to limit data
  • slice_width (integer) – width of the image slice, default 10 pixels
Returns:

top, bottom, left, right – image slices oriented in the same direction

Return type:

1D np.arrays shape (rmin:rmax, 1)

abel.tools.center.find_image_center_by_slice(IM, slice_width=10, radial_range=(0, -1), axis=(0, 1), **kwargs)[source]

Center image by comparing opposite side, vertical (axis=0) and/or horizontal slice (axis=1) profiles. To center along both axis, use axis=(0,1).

Parameters:
  • IM (2D np.array) – The image data.
  • slice_width (integer) – Sum together this number of rows (cols) to improve signal, default 10.
  • radial_range (tuple) – (rmin,rmax): radial range [rmin:rmax] for slice profile comparison.
  • axis (integer or tuple) – Center with along axis=0 (vertical), or axis=1 (horizontal), or axis=(0,1) (both vertical and horizontal).
Returns:

(vertical_shift, horizontal_shift) – (axis=0 shift, axis=1 shift)

Return type:

tuple of floats

abel.tools.circularize module

abel.tools.circularize.circularize_image(IM, method='lsq', center=None, radial_range=None, dr=0.5, dt=0.5, smooth=0, ref_angle=None, inverse=False, return_correction=False)[source]

Corrects image distortion on the basis that the structure should be circular.

This is a simplified radial scaling version of the algorithm described in J. R. Gascooke and S. T. Gibson and W. D. Lawrance: ‘A “circularisation” method to repair deformations and determine the centre of velocity map images’ J. Chem. Phys. 147, 013924 (2017).

This function is especially useful for correcting the image obtained with a velocity-map-imaging spectrometer, in the case where there is distortion of the Newton Sphere (ring) structure due to an imperfect electrostatic lens or stray electromagnetic fields. The correction allows the highest-resolution 1D photoelectron distribution to be extracted.

The algorithm splits the image into “slices” at many different angles (set by dt) and compares the radial intensity profile of adjacent slices. A scaling factor is found which aligns each slice profile with the previous slice. The image is then corrected using a spline function that smoothly connects the discrete scaling factors as a continuous function of angle.

This circularization algorithm should only be applied to a well-centered image, otherwise use the center keyword (described below) to center it.

Parameters:
  • IM (numpy 2D array) – Image to be circularized.
  • method (str) –

    Method used to determine the radial correction factor to align slice profiles:

    argmax - compare intensity-profile.argmax() of each radial slice.
    This method is quick and reliable, but it assumes that the radial intensity profile has an obvious maximum. The positioning is limited to the nearest pixel.
    lsq - minimize the difference between a slice intensity-profile
    with its adjacent slice. This method is slower and may fail to converge, but it may be applied to images with any (circular) structure. It aligns the slices with sub-pixel precision.
  • center (str, float tuple, or None) – Pre-center image using abel.tools.center.center_image(). center may be: com, convolution, gaussian, image_center, slice, or a float tuple center \((y, x)\).
  • radial_range (tuple, or None) – Limit slice comparison to the radial range tuple (rmin, rmax), in pixels, from the image center. Use to determine the distortion correction associated with particular peaks. It is recommended to select a region of your image where the signal-to-noise is highest, with sharp persistant (in angle) features.
  • dr (float) –

    Radial grid size for the polar coordinate image, default = 0.5 pixel. This is passed to abel.tools.polar.reproject_image_into_polar().

    Small values may improve the distortion correction, which is often of sub-pixel dimensions, at the cost of reduced signal to noise for the slice intensity profile. As a general rule, dr should be significantly smaller than the radial “feature size” in the image.

  • dt (float) –

    Angular grid size. This sets the number of radial slices, given by \(2\pi/dt\). Default = 0.1, ~ 63 slices. More slices, using smaller dt, may provide a more detailed angular variation of the correction, at the cost of greater signal to noise in the correction function.

    Also passed to abel.tools.polar.reproject_image_into_polar()

  • smooth (float) –

    This value is passed to the scipy.interpolate.UnivariateSpline() function and controls how smooth the spline interpolation is. A value of zero corresponds to a spline that runs through all of the points, and higher values correspond to a smoother spline function.

    It is important to examine the relative peak position (scaling factor) data and how well it is represented by the spline function. Use the option return_correction=True to examine this data. Typically, smooth may remain zero, noisy data may require some smoothing.

  • ref_angle (None or float) –

    Reference angle for which radial coordinate is unchanged. Angle varies between \(-\pi\) to \(\pi\), with zero angle vertical.

    None uses numpy.mean(radial scale factors)(), which attempts to maintain the same average radial scaling. This approximation is likely valid, unless you know for certain that a specific angle of your image corresponds to an undistorted image.

  • inverse (bool) –

    Apply an inverse Abel transform the polar coordinate image, to remove the background intensity. This may improve the signal to noise, allowing the weaker intensity featured to be followed in angle.

    Note that this step is only for the purposes of allowing the algorithm to better follow peaks in the image. It does not affect the final image that is returned, except for (hopefully) slightly improving the precision of the distortion correction.

  • return_correction (bool) – Additional outputs, as describe below.
Returns:

  • IMcirc (numpy 2D array, same size as input) – Circularized version of the input image.

    The following values are returned if return_correction=True:

  • angles (numpy 1D array) – Mid-point angle (radians) of each image slice.

  • radial_correction (numpy 1D array) – Radial correction scale factor at each angular slice.

  • radial_correction_function (numpy function that accepts numpy.array) – Function that may be used to evaluate the radial correction at any angle.

abel.tools.circularize.circularize(IM, radial_correction_function, ref_angle=None)[source]

Remap image from its distorted grid to the true cartesian grid.

Parameters:
  • IM (numpy 2D array) – Original image
  • radial_correction_function (funct) – A function returning the radial correction for a given angle. It should accept a numpy 1D array of angles.
abel.tools.circularize.correction(polarIMTrans, angles, radial, method)[source]
Determines a radial correction factors that align an angular slice
radial intensity profile with its adjacent (previous) slice profile.
Parameters:
  • polarIMTrans (numpy 2D array) – Polar coordinate image, transposed \((\theta, r)\) so that each row is a single angle.
  • angles (numpy 1D array) – Angle coordinates for one row of polarIMTrans.
  • radial (numpy 1D array) – Radial coordinates for one column of polarIMTrans.
  • method (str) –

    “argmax”: radial correction factor from position of maximum intensity.

    ”lsq” : least-squares determine a radial correction factor that will align a radial intensity profile with the previous, adjacent slice.

abel.tools.math module

abel.tools.math.gradient(f, x=None, dx=1, axis=-1)[source]

Return the gradient of 1 or 2-dimensional array. The gradient is computed using central differences in the interior and first differences at the boundaries.

Irregular sampling is supported (it isn’t supported by np.gradient)

Parameters:
  • f (1d or 2d numpy array) – Input array.
  • x (array_like, optional) – Points where the function f is evaluated. It must be of the same length as f.shape[axis]. If None, regular sampling is assumed (see dx)
  • dx (float, optional) – If x is None, spacing given by dx is assumed. Default is 1.
  • axis (int, optional) – The axis along which the difference is taken.
Returns:

out – Returns the gradient along the given axis.

Return type:

array_like

Notes

To-Do: implement smooth noise-robust differentiators for use on experimental data. http://www.holoborodko.com/pavel/numerical-methods/numerical-derivative/smooth-low-noise-differentiators/

abel.tools.math.gaussian(x, a, mu, sigma, c)[source]

Gaussian function

\(f(x)=a e^{-(x - \mu)^2 / (2 \sigma^2)} + c\)

ref: https://en.wikipedia.org/wiki/Gaussian_function

Parameters:
  • x (1D np.array) – coordinate
  • a (float) – the height of the curve’s peak
  • mu (float) – the position of the center of the peak
  • sigma (float) – the standard deviation, sometimes called the Gaussian RMS width
  • c (float) – non-zero background
Returns:

out – the Gaussian profile

Return type:

1D np.array

abel.tools.math.guss_gaussian(x)[source]

Find a set of better starting parameters for Gaussian function fitting

Parameters:x (1D np.array) – 1D profile of your data
Returns:out – estimated value of (a, mu, sigma, c)
Return type:tuple of float
abel.tools.math.fit_gaussian(x)[source]

Fit a Gaussian function to x and return its parameters

Parameters:x (1D np.array) – 1D profile of your data
Returns:out – (a, mu, sigma, c)
Return type:tuple of float

abel.tools.polar module

abel.tools.polar.reproject_image_into_polar(data, origin=None, Jacobian=False, dr=1, dt=None)[source]

Reprojects a 2D numpy array (data) into a polar coordinate system. “origin” is a tuple of (x0, y0) relative to the bottom-left image corner, and defaults to the center of the image.

Parameters:
  • data (2D np.array) –
  • origin (tuple) – The coordinate of the image center, relative to bottom-left
  • Jacobian (boolean) – Include r intensity scaling in the coordinate transform. This should be included to account for the changing pixel size that occurs during the transform.
  • dr (float) – Radial coordinate spacing for the grid interpolation tests show that there is not much point in going below 0.5
  • dt (float) – Angular coordinate spacing (in radians) if dt=None, dt will be set such that the number of theta values is equal to the maximum value between the height or the width of the image.
Returns:

  • output (2D np.array) – The polar image (r, theta)
  • r_grid (2D np.array) – meshgrid of radial coordinates
  • theta_grid (2D np.array) – meshgrid of theta coordinates

Notes

Adapted from: http://stackoverflow.com/questions/3798333/image-information-along-a-polar-coordinate-system

abel.tools.polar.index_coords(data, origin=None)[source]

Creates x & y coords for the indicies in a numpy array

Parameters:
  • data (numpy array) – 2D data
  • origin ((x,y) tuple) – defaults to the center of the image. Specify origin=(0,0) to set the origin to the bottom-left corner of the image.
Returns:

x, y

Return type:

arrays

abel.tools.polar.cart2polar(x, y)[source]

Transform Cartesian coordinates to polar

Parameters:y (x,) – Cartesian coordinates
Returns:r, theta – Polar coordinates
Return type:floats or arrays
abel.tools.polar.polar2cart(r, theta)[source]

Transform polar coordinates to Cartesian

Parameters:theta (r,) – Polar coordinates
Returns:x, y – Cartesian coordinates
Return type:floats or arrays

abel.tools.transform_pairs module

Analytical function Abel transform pairs

profiles 1-7, table 1 of: G. C.-Y Chan and G. M. Hieftje Spectrochimica Acta B 61, 31-41 (2006)

Note: profile4 does not produce a correct Abel transform pair due to typographical errors in the publications

profile 8, curve B in table 2 of: Hansen and Law J. Opt. Soc. Am. A 2 510-520 (1985)

Note: the transform pair functions are more conveniently accessed via the class:

   func = abel.tools.analytical.TransformPair(n, profile=nprofile)

which sets the radial range r and provides attributes:
``.func`` (source), ``.abel`` (projection), ``.r`` (radial range),
``.dr`` (step), ``.label`` (the profile name)
r : floats or numpy 1D array of floats
value or grid to evaluate the function pair: 0 < r < 1
source, projection : tuple of 1D numpy arrays of shape r
source function profile (inverse Abel transform of projection), projection functon profile (forward Abel transform of source)
abel.tools.transform_pairs.a(n, r)[source]

coefficient

\[a_n = \sqrt{n^2 - r^2}\]
abel.tools.transform_pairs.profile1(r)[source]

profile1: Cremers and Birkebak App. Opt. 5, 1057-1064 (1966) Eq(13)

\[ \begin{align}\begin{aligned}\epsilon(r) &= 0.75 + 12r^2 - 32r^3 & 0 \le r \le 0.25\\\epsilon(r) &= \frac{16}{27}(1 + 6r - 15r^2 + 8r^3) & 0.25 \lt r \le 1\\I(r) &= \frac{1}{108}(128a_1 +a_{0.25}) + \frac{2}{27}r^2 (283a_{0.25} - 112a_1) +\\& \,\,\,\, \frac{8}{9}r^2\left[4(1+r^2)\ln\frac{1+a_1}{r} - (4+31r^2)\ln\frac{0.25+a_{0.25}}{r}\right] & 0 \le r \le 0.25\\I(r) &= \frac{32}{27}\left[a_1 - 7a_1 r + 3r^2(1+r^2) \ln\frac{1+a_1}{r}\right] & 0.25 \lt r \le 1\end{aligned}\end{align} \]
               profile1
      source                projection
┼+1.3                  ┼+1.3
│                      o   o
│     x                │     o
│   x   x              │       o
│ x                    │
x         x            │         o
│           x          │
│                      │           o
│             x        │
│                      │             o
┼+0─────────────x──┼   ┼+0─────────────o──┼
0          r      +1   0          r      +1
abel.tools.transform_pairs.profile2(r)[source]

profile2: Cremers and Birkebak App. Opt. 5, 1057-1064 (1966) Eq(13)

\[ \begin{align}\begin{aligned}\epsilon(r) &= 1 - 3r^2 + 2r^3 & 0 \le r \le 1\\I(r) &= a_1\left(1-\frac{5}{2}r^2\right) + \frac{3}{2}r^4\ln\frac{1+a_1}{r} & 0 \le r \le 1\end{aligned}\end{align} \]
               profile2
      source                projection
┼+1.1                  ┼+1.1
x x                    o o
│   x                  │   o
│     x                │     o
│       x              │
│                      │       o
│         x            │
│                      │         o
│           x          │           o
│             x        │
┼+0─────────────x──┼   ┼+0───────────o────┼
0          r      +1   0          r      +1
abel.tools.transform_pairs.profile3(r)[source]

profile3: Cremers and Birkebak App. Opt. 5, 1057-1064 (1966) Eq(13)

\[ \begin{align}\begin{aligned}\epsilon(r) &= 1-2r^2 & 0 \le r \le 0.5\\\epsilon(r) &= 2(1-r^2)^2 & 0.5 \lt r \le 1\\I(r) &= \frac{4a_1}{3}(1+2r^2)-\frac{2 a_{0.5}}{3}(1+8r^2) - 4r^2\ln\frac{1-a_1}{0.5+a_{0.5}} & 0 \le r \le 0.5\\I(r) &= \frac{4a_1}{3}(1+2r^2)-4r^2\ln\frac{1-a_1}{r} & 0.5 \lt r \le 1\end{aligned}\end{align} \]
               profile3
      source                projection
┼+1.1                  ┼+1.1
x x x                  o o
│                      │   o
│     x                │     o
│       x              │
│                      │       o
│         x            │
│                      │         o
│           x          │
│                      │           o
┼+0───────────x────┼   ┼+0───────────o────┼
0          r      +1   0          r      +1
abel.tools.transform_pairs.profile4(r)[source]
profile4: `Alvarez, Rodero, Quintero Spectochim. Acta B 57,

1665-1680 (2002) <https://doi.org/10.1016/S0584-8547(02)00087-3>`_

WARNING: projection function pair incorrect due to typo errors in Table 1.

\[ \begin{align}\begin{aligned}\epsilon(r) &= 0.1 + 5.5r^2 - 5.25r^3 & 0 \le r \le 0.7\\\epsilon(r) &= -40.74 + 155.56r - 188.89r^2 + 74.07r^3 & 0.7 \lt r \le1\\I(r) &= 22.68862a_{0.7} - 14.811667a_1 + (217.557a_{0.7} - 193.30083a_1)r^2 +\\ & \,\,\, 155.56r^2\ln\frac{1 + a_1}{0.7 + a_{0.7}} + r^4\left(55.5525\ln\frac{1 + a_1}{r} - 59.49\ln\frac{0.7 + a_{0.7}}{r}\right) & 0 \le r \le 0.7\\I(r) &= -14.811667a_1 - 193.30083a_1 r^2 + r^2(155.56 + 55.5525r^2) \ln\frac{1 + a_1}{r} & 0.7 \lt r \le 1\end{aligned}\end{align} \]
               profile4     (incorrect)
      source                projection
┼+2.2                  ┼+2.2       o
│                      │         o   o
│                      │       o       o
│                      │     o
│                      │
│                      │ o o             o
│           x x        o
│         x     x      │
│       x              │
│     x                │
┼+0─x─────────────x┼   ┼+0────────────────┼
0          r      +1   0          r      +1
abel.tools.transform_pairs.profile5(r)[source]

profile5: Buie et al. J. Quant. Spectrosc. Radiat. Transfer 55, 231-243 (1996)

\[ \begin{align}\begin{aligned}\epsilon(r) &= 1 & 0 \le r \le 1\\I(r) &= 2a_1 & 0 \le r \le 1\end{aligned}\end{align} \]
               profile5
      source                projection
┼+2.1                  ┼+2.1
│                      │     o o
│                      │         o o
│                      │             o
│                      │
│                      │               o
x x x x x x x x x x    │
│                      │                 o
│                      │
│                      │
┼+0────────────────┼   ┼+0────────────────┼
0          r      +1   0          r      +1
abel.tools.transform_pairs.profile6(r)[source]

profile6: Buie et al. J. Quant. Spectrosc. Radiat. Transfer 55, 231-243 (1996)

\[ \begin{align}\begin{aligned}\epsilon(r) &= (1-r^2)^{-\frac{3}{2}} \exp\left[1.1^2\left( 1 - \frac{1}{1-r^2}\right)\right] & 0 \le r \le 1\\I(r) &= \frac{\sqrt{\pi}}{1.1a_1} \exp\left[1.1^2\left( 1 - \frac{1}{1-r^2}\right)\right] & 0 \le r \le 1\end{aligned}\end{align} \]
               profile6
      source                projection
┼+1.8                  ┼+1.8
│                      o o o
│                      │     o o
│                      │         o
│                      │
x x x x x x x          │           o
│             x        │
│                      │             o
│               x      │
│                      │               o
┼+0────────────────┼   ┼+0────────────────┼
0          r      +1   0          r      +1
abel.tools.transform_pairs.profile7(r)[source]

profile7: Buie et al. J. Quant. Spectrosc. Radiat. Transfer 55, 231-243 (1996)

\[ \begin{align}\begin{aligned}\epsilon(r) &= \frac{1}{2}(1+10r^2-23r^4+12r^6) & 0 \le r \le 1\\I(r) &= \frac{8}{105}a_1(19 + 34r^2 - 125r^4 + 72r^6) & 0 \le r \le 1\end{aligned}\end{align} \]
               profile7
      source                projection
┼+1.7                  ┼+1.7
│                      o o o o o
│                      │         o
│                      │
│       x x x          │           o
│     x       x        │
│                      │             o
│   x                  │
x x             x      │
│                      │               o
┼+0───────────────x┼   ┼+0────────────────┼
0          r      +1   0          r      +1
abel.tools.transform_pairs.profile8(r)[source]

profile8: Curve B table 2 of Hansen and Law J. Opt. Soc. Am. A 2 510-520 (1985)

\[ \begin{align}\begin{aligned}\epsilon(r) &= (1-r^2)^{-\frac{3}{2}} \exp\left[\frac{(1.1r)^2}{r^2-1}\right]\\I(r) &= \frac{\pi^\frac{1}{2}}{1.1}(1-r^2)^{-\frac{1}{2}} \exp\left[\frac{(1.1r)^2}{r^2-1}\right]\end{aligned}\end{align} \]
               profile8
      source                projection
┼+1.8                  ┼+1.8
│                      o o o
│                      │     o o
│                      │         o
│                      │
x x x x x x x          │           o
│             x        │
│                      │             o
│               x      │
│                      │               o
┼+0────────────────┼   ┼+0────────────────┼
0          r      +1   0          r      +1

abel.tools.symmetry module

abel.tools.symmetry.get_image_quadrants(IM, reorient=True, symmetry_axis=None, use_quadrants=(True, True, True, True), symmetrize_method=u'average')[source]

Given an image (m,n) return its 4 quadrants Q0, Q1, Q2, Q3 as defined below.

Parameters:
  • IM (2D np.array) – Image data shape (rows, cols)
  • reorient (boolean) – Reorient quadrants to match the orientation of Q0 (top-right)
  • symmetry_axis (int or tuple) – can have values of None, 0, 1, or (0,1) and specifies no symmetry, vertical symmetry axis, horizontal symmetry axis, and both vertical and horizontal symmetry axes. Quadrants are added. See Note.
  • use_quadrants (boolean tuple) – Include quadrant (Q0, Q1, Q2, Q3) in the symmetry combination(s) and final image
  • symmetrize_method (str) –

    Method used for symmetrizing the image.

    average: Simply average the quadrants. fourier: Axial symmetry implies that the Fourier components of the 2-D

    projection should be real. Removing the imaginary components in reciprocal space leaves a symmetric projection. ref: Overstreet, K., et al. “Multiple scattering and the density
    distribution of a Cs MOT.” Optics express 13.24 (2005): 9672-9682. http://dx.doi.org/10.1364/OPEX.13.009672
Returns:

Q0, Q1, Q2, Q3 – shape: (rows//2+rows%2, cols//2+cols%2) all oriented in the same direction as Q0 if reorient=True

Return type:

tuple of 2D np.arrays

Notes

The symmetry_axis keyword averages quadrants like this:

 +--------+--------+
 | Q1   * | *   Q0 |
 |   *    |    *   |
 |  *     |     *  |               cQ1 | cQ0
 +--------o--------+ --(output) -> ----o----
 |  *     |     *  |               cQ2 | cQ3
 |   *    |    *   |
 | Q2  *  | *   Q3 |          cQi == combined quadrants
 +--------+--------+

symmetry_axis = None - individual quadrants
symmetry_axis = 0 (vertical) - average Q0+Q1, and Q2+Q3
symmetry_axis = 1 (horizontal) - average Q1+Q2, and Q0+Q3
symmetry_axis = (0, 1) (both) - combine and average all 4 quadrants

The end results look like this:

(0) symmetry_axis = None

    returned image   Q1 | Q0
                    ----o----
                     Q2 | Q3

(1) symmetry_axis = 0

    Combine:  Q01 = Q0 + Q1, Q23 = Q2 + Q3
    returned image    Q01 | Q01
                     -----o-----
                      Q23 | Q23

(2) symmetry_axis = 1

    Combine: Q12 = Q1 + Q2, Q03 = Q0 + Q3
    returned image   Q12 | Q03
                    -----o-----
                     Q12 | Q03

(3) symmetry_axis = (0, 1)

    Combine all quadrants: Q = Q0 + Q1 + Q2 + Q3
    returned image   Q | Q
                    ---o---  all quadrants equivalent
                     Q | Q
abel.tools.symmetry.put_image_quadrants(Q, original_image_shape, symmetry_axis=None)[source]

Reassemble image from 4 quadrants Q = (Q0, Q1, Q2, Q3) The reverse process to get_image_quadrants(reorient=True)

Note: the quadrants should all be oriented as Q0, the upper right quadrant

Parameters:
  • Q (tuple of np.array (Q0, Q1, Q2, Q3)) –

    Image quadrants all oriented as Q0 shape (rows//2+rows%2, cols//2+cols%2)

    +--------+--------+
    | Q1   * | *   Q0 |
    |   *    |    *   |
    |  *     |     *  |
    +--------o--------+
    |  *     |     *  |
    |   *    |    *   |
    | Q2  *  | *   Q3 |
    +--------+--------+
    
  • original_image_shape (tuple) –

    (rows, cols)

    reverses the padding added by get_image_quadrants() for odd-axis sizes

    odd row trims 1 row from Q1, Q0

    odd column trims 1 column from Q1, Q2

  • symmetry_axis (int or tuple) –

    impose image symmetry

    symmetry_axis = 0 (vertical)   - Q0 == Q1 and Q3 == Q2 symmetry_axis = 1 (horizontal) - Q2 == Q1 and Q3 == Q0
Returns:

IM

Reassembled image of shape (rows, cols):

 symmetry_axis =

  None             0              1           (0,1)

 Q1 | Q0        Q1 | Q1        Q1 | Q0       Q1 | Q1
----o----  or  ----o----  or  ----o----  or ----o----
 Q2 | Q3        Q2 | Q2        Q1 | Q0       Q1 | Q1

Return type:

np.array

abel.tools.vmi module

abel.tools.vmi.angular_integration(IM, origin=None, Jacobian=True, dr=1, dt=None)[source]

Angular integration of the image.

Returns the one-dimensional intensity profile as a function of the radial coordinate.

Note: the use of Jacobian=True applies the correct Jacobian for the integration of a 3D object in spherical coordinates.

Parameters:
  • IM (2D numpy.array) – The data image.
  • origin (tuple) – Image center coordinate relative to bottom-left corner defaults to rows//2, cols//2.
  • Jacobian (boolean) – Include \(r\sin\theta\) in the angular sum (integration). Also, Jacobian=True is passed to abel.tools.polar.reproject_image_into_polar(), which includes another value of r, thus providing the appropriate total Jacobian of \(r^2\sin\theta\).
  • dr (float) – Radial coordinate grid spacing, in pixels (default 1). dr=0.5 may reduce pixel granularity of the speed profile.
  • dt (float) – Theta coordinate grid spacing in radians. if dt=None, dt will be set such that the number of theta values is equal to the height of the image (which should typically ensure good sampling.)
Returns:

  • r (1D numpy.array) – radial coordinates
  • speeds (1D numpy.array) – Integrated intensity array (vs radius).

abel.tools.vmi.average_radial_intensity(IM, **kwargs)[source]

Calculate the average radial intensity of the image, averaged over all angles. This differs form abel.tools.vmi.angular_integration() only in that it returns the average intensity, and not the integrated intensity of a 3D image. It is equivalent to calling abel.tools.vmi.angular_integration() with Jacobian=True and then dividing the result by 2*pi.

Parameters:
Returns:

  • r (1D numpy.array) – radial coordinates
  • intensity (1D numpy.array) – one-dimensional intensity profile as a function of the radial coordinate.

abel.tools.vmi.radial_integration(IM, radial_ranges=None)[source]

Intensity variation in the angular coordinate.

This function is the \(\theta\)-coordinate complement to abel.tools.vmi.angular_integration()

Evaluates intensity vs angle for defined radial ranges. Determines the anisotropy parameter for each radial range.

See examples/example_PAD.py

Parameters:
  • IM (2D numpy.array) – Image data
  • radial_ranges (list of tuple ranges or int step) –
    tuple integration ranges
    [(r0, r1), (r2, r3), ...] evaluates the intensity vs angle for the radial ranges r0_r1, r2_r3, etc.

    int - the whole radial range (0, step), (step, 2*step), ..

Returns:

  • Beta (array of tuples) – (beta0, error_beta_fit0), (beta1, error_beta_fit1), … corresponding to the radial ranges
  • Amplitude (array of tuples) – (amp0, error_amp_fit0), (amp1, error_amp_fit1), … corresponding to the radial ranges
  • Rmidpt (numpy float 1d array) – radial-mid point of each radial range
  • Intensity_vs_theta (2D numpy.array) – Intensity vs angle distribution for each selected radial range.
  • theta (1D numpy.array) – Angle coordinates, referenced to vertical direction.

abel.tools.vmi.anisotropy_parameter(theta, intensity, theta_ranges=None)[source]

Evaluate anisotropy parameter \(\beta\), for \(I\) vs \(\theta\) data.

\[I = \frac{\sigma_\text{total}}{4\pi} [ 1 + \beta P_2(\cos\theta) ]\]

where \(P_2(x)=\frac{3x^2-1}{2}\) is a 2nd order Legendre polynomial.

Cooper and Zare “Angular distribution of photoelectrons” J Chem Phys 48, 942-943 (1968)

Parameters:
  • theta (1D numpy array) – Angle coordinates, referenced to the vertical direction.
  • intensity (1D numpy array) – Intensity variation with angle
  • theta_ranges (list of tuples) – Angular ranges over which to fit [(theta1, theta2), (theta3, theta4)]. Allows data to be excluded from fit, default include all data
Returns:

  • beta (tuple of floats) – (anisotropy parameter, fit error)
  • amplitude (tuple of floats) – (amplitude of signal, fit error)

abel.tools.vmi.toPES(radial, intensity, energy_cal_factor, per_energy_scaling=True, photon_energy=None, Vrep=None, zoom=1)[source]

Convert speed radial coordinate into electron kinetic or electron binding energy. Return the photoelectron spectrum (PES).

This calculation uses a single scaling factor energy_cal_factor to convert the radial pixel coordinate into electron kinetic energy.

Additional experimental parameters: photon_energy will give the energy scale as electron binding energy, in the same energy units, while Vrep, the VMI lens repeller voltage (volts), provides for a voltage independent scaling factor. i.e. energy_cal_factor should remain approximately constant.

The energy_cal_factor is readily determined by comparing the generated energy scale with published spectra. e.g. for O- photodetachment, the strongest fine-structure transition occurs at the electron affinity \(EA = 11,784.676(7)\) cm \(^{-1}\). Values for the ANU experiment are given below, see also examples/example_hansenlaw.py.

Parameters:
  • radial (numpy 1D array) – radial coordinates.
  • intensity (numpy 1D array) – intensity values, at the radial array.
  • energy_cal_factor (float) – energy calibration factor that will convert radius squared into energy. The units affect the units of the output. e.g. inputs in eV/pixel2, will give output energy units in eV. A value of \(1.148427\times 10^{-5}\) cm\(^{-1}/\)pixel2 applies for “examples/data/O-ANU1024.txt” (with Vrep = -98 volts).
  • per_energy_scaling (bool) –
    sets the intensity Jacobian.
    If True, the returned intensities correspond to an “intensity per eV” or “intensity per cm-1 “. If False, the returned intensities correspond to an “intensity per pixel”.

    Optional:

  • photon_energy (None or float) – measurement photon energy. The output energy scale is then set to electron-binding-energy in units of energy_cal_factor. The conversion from wavelength (nm) to photon_energy in (cm:sup:-1) is \(10^{7}/\lambda\) (nm) e.g. 1.0e7/812.51 for “examples/data/O-ANU1024.txt”.
  • Vrep (None or float) – repeller voltage. Convenience parameter to allow the energy_cal_factor to remain constant, for different VMI lens repeller voltages. Defaults to None, in which case no extra scaling is applied. e.g. -98 volts, for “examples/data/O-ANU1024.txt”.
  • zoom (float) – additional scaling factor if the input experimental image has been zoomed. Default 1.
Returns:

  • eKBE (numpy 1d-array of floats) – energy scale for the photoelectron spectrum in units of energy_cal_factor. Note that the data is no-longer on a uniform grid.
  • PES (numpy 1d-array of floats) – the photoelectron spectrum, scaled according to the per_energy_scaling input parameter.

abel.benchmark module

class abel.benchmark.AbelTiming(n=[301, 501], select=[u'all'], n_max_bs=700, n_max_slow=700, transform_repeat=1)[source]

Bases: object

abel.benchmark.is_symmetric(arr, i_sym=True, j_sym=True)[source]

Takes in an array of shape (n, m) and check if it is symmetric

Parameters:
  • arr (1D or 2D array) –
  • i_sym (array) – symmetric with respect to the 1st axis
  • j_sym (array) – symmetric with respect to the 2nd axis
Returns:

  • a binary array with the symmetry condition for the corresponding quadrants.
  • The globa
  • Note (if both i_sym=True and i_sym=True, the input array is checked)
  • for polar symmetry.
  • See https (//github.com/PyAbel/PyAbel/issues/34#issuecomment-160344809)
  • for the defintion of a center of the image.

abel.benchmark.absolute_ratio_benchmark(analytical, recon, kind=u'inverse')[source]

Check the absolute ratio between an analytical function and the result of a inv. Abel reconstruction.

Parameters:
  • analytical (one of the classes from analytical, initialized) –
  • recon (1D ndarray) – a reconstruction (i.e. inverse abel) given by some PyAbel implementation

abel.tests module

abel.tests.run.run(coverage=False)[source]

This runs the complete set of PyAbel tests.

abel.tests.run.run_cli(coverage=False)[source]

This also runs the complete set of PyAbel tests. But uses sys.exit(status) instead of simply returning the status.