abel package

abel.transform module

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

This is the main function of PyAbel, providing both forward and inverse abel transforms for full images. In addition, this function can perform image centering and symmetrization.

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.
    three_point
    the three-point transform of Dasch and co-workers
  • 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.
    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.
  • verbose (boolean) – True/False to determine if non-critical output should be printed.
  • symmetry_axis (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.
  • 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.)
  • 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.

Note

Quadrant averaging

The quadrants can be 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-----
                   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
    
Returns:results – The transform function returns a dictionary of results depending on the options selected
results['transform']
(always returned) is the 2D forward/reverse Abel transform
results['radial_intensity']
is not currently implemented
results['residual']
is not currently implemented
Return type:dict

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.
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.

abel.basex module

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.
  • 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.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')

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

Note: only the 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.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.hansenlaw module

abel.hansenlaw.hansenlaw_transform(IM, dr=1, direction=u'inverse')[source]

Forward/Inverse Abel transformation using the algorithm of Hansen and Law J. Opt. Soc. Am. A 2, 510-520 (1985) equation 2a:

\[f(r) = \frac{1}{\pi} \int_{r}^{\inf} \frac{g^\prime(R)}{\sqrt{R^2-r^2}} dR,\]

where

\(f(r)\) is the reconstructed image (source) function, \(g'(R)\) is the derivative of the projection (measured) function

Evaluation follows Eqs. (15 or 17), using (16a), (16b), and (16c or 18) of the Hansen and Law paper. For the full image transform, use abel.transform.

For the inverse Abel transform of image g:

f = abel.transform(g, direction="inverse", method="hansenlaw")["transform"]

For the forward Abel transform of image f:

g = abel.transform(r, direction="forward", method="hansenlaw")["transform"]

This function performs the Hansen-Law transform on only one “right-side” image, typically one quadrant of the full image:

Qtrans = abel.hansenlaw.hansenlaw_transform(Q, direction="inverse")

Recursion method proceeds from the outer edge of the image toward the image centre (origin). i.e. when n=N-1, R=Rmax, and when n=0, R=0. This fits well with processing the image one quadrant (chosen orientation to be rightside-top), or one right-half image at a time.

Parameters:
  • IM (2D np.array) – One quadrant (or half) of the image oriented top-right.
  • dr (float) – Sampling size (=1 for pixel images), used for Jacobian scaling
  • direction (string ('forward' or 'inverse')) – forward or inverse Abel transform
Returns:

AIM – forward/inverse Abel transform image

Return type:

2D numpy array

Note

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

.         +--------      +--------+
.         |      *       | *      |
.         |   *          |    *   |  <---------- IM
.         |  *           |     *  |
.         +--------      o--------+
.         |  *           |     *  |
.         |   *          |    *   |
.         |     *        | *      |
.         +--------      +--------+

Image centre o should be within a pixel (i.e. an odd number of columns) Use abel.tools.center.center_image(IM, method='com', odd_size=True)

abel.three_point module

abel.three_point.OP0(i, j)[source]

This is the I_ij(0) function in Dasch 1992, pg 1147, Eq (7)

abel.three_point.OP1(i, j)[source]

This is the I_ij(1) function in Dasch 1992, pg 1147, Eq (7)

abel.three_point.OP_D(i, j)[source]

Calculates three-point abel inversion operator D_ij, following Eq (6) in Dasch 1992 (Applied Optics). The original reference contains several typos. One correction is done in function OP1 following Karl Martin’s PhD thesis See here: https://www.lib.utexas.edu/etd/d/2002/martinkm07836/martinkm07836.pdf

abel.three_point.get_bs_three_point_cached(col, basis_dir='.', verbose=False)[source]

Internal function.

Gets thre_point operator matrix corresponding to specified image size, using the disk as a cache. (i.e., load from disk if they exist, if not, calculate them and save a copy on disk)

Parameters:
  • col (integer) – Width of image to be inverted using three_point method. Three_point operator matrix will be of size (col x col)
  • basis_dir (string) – path to the directory for saving / loading the three_point operator matrix. If None, the operator matrix will not be saved to disk.
  • verbose (True/False) – Set to True to see more output for debugging
abel.three_point.three_point(data, center, dr=1.0, basis_dir='./', verbose=False, direction='inverse')[source]

This function splits the image into two halves, sends each half to three_point_transform(), stitches the output back together, and returns the full transform to the user.

Parameters:
  • data (NxM numpy array) – The raw data is presumed to be symmetric about the vertical axis.
  • center (integer or tuple (x,y)) – The location of the symmetry axis (center column index of the image) or the center of the image in (x,y) format.
  • dr (float) – Size of one pixel in the radial direction
  • basis_dir (string) – path to the directory for saving / loading the three_point operator matrix. If None, the operator matrix will not be saved to disk.
  • verbose (True/False) – Set to True to see more output for debugging
  • direction (str) – The type of Abel transform to be performed. Currently only accepts value ‘inverse’
Returns:

inv_IM – Abel inversion of IM - a rows x cols array

Return type:

numpy array

abel.three_point.three_point_core_transform(IM, D)[source]

This is the internal function that performs the actual three_point transform. It requires that the matrix of basis set coefficients, D, be passed as a parameter.

Parameters:
  • IM (numpy array) – Raw data - a rows x cols array
  • D (basis sets) – generated by the get_bs_three_point_cached() function
Returns:

inv_IM – Abel inversion of IM - a rows x cols array

Return type:

numpy array

abel.three_point.three_point_transform(IM, basis_dir='.', direction='inverse', verbose=False)[source]

Inverse Abel transformation using the algorithm of: Dasch, Applied Optics, Vol. 31, No. 8, 1146-1152 (1992).

Parameters:IM (numpy array) – Raw data - a rows x cols array
Returns:inv_IM – Abel inversion of IM - a rows x cols array
Return type:numpy array

abel.direct module

abel.direct.simpson_rule_wrong(f, x=None, dx=None, axis=1, **args)[source]

This function computes the Simpson rule https://en.wikipedia.org/wiki/Simpson%27s_rule#Python both in the cases of odd and even intervals which is not technically valid

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

This algorithm does a direct computation of the Abel transform

  • integration near the singular value is done analytically
  • integration further from the singular value with the Simpson rule.
Parameters:
  • fr (1d or 2d numpy array) – input array to which direct/inversed 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)
  • f (1D ndarray) – the spatial mesh (optional)
  • derivative (callable) – a function that can return the derivative of the fr array with respect to r (only used in the inverse Abel transform).
  • inverse (boolean) – If True inverse Abel transform is applied, otherwise use a forward Abel transform.
  • correction (boolean) – if False integration is performed with the Simpson rule, the pixel where the weighting function has a singular value is ignored if True in addition to the integration with the Simpson rule, integration near the singular value is done analytically, assuming a piecewise linear data.
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.

abel.direct.reflect_array(x, axis=1, kind=u'even')[source]

Make a symmetrically reflected array with respect to the given axis

abel.tools.analytical module

abel.tools.analytical.sample_image(n=361, name='dribinski', sigma=2, temperature=200)[source]

Sample images, made up of Gaussian functions

Parameters:
  • n (integer) – image size n rows x n cols
  • name (str) – one of “dribinski” or “Ominus”
  • sigma (float) – Gaussian 1/e width (pixels)
  • temperature (float) – for ‘Ominus’ only anion levels have Boltzmann population weight (2J+1) exp(-177.1 h c 100/k/temperature)
Returns:

IM – image

Return type:

2D np.array

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

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

abel.tools.analytical.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

This function is mostly used for unit testing the inverse Abel transform

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

abel.tools.analytical.sym_abel_step_1d(r, A0, r0, r1)[source]

Produces a symmetrical analytical transform of a 1d step

abel.tools.center module

abel.tools.center.find_center(IM, method=u'image_center', verbose=False, **kwargs)[source]
Parameters:
  • IM (2D np.array) – image data
  • method (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
    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.
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, 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
    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.
  • 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.
Returns:

out – Centered image

Return type:

2D np.array

abel.tools.center.set_center(data, center, crop=u'maintain_size', 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.
  • 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_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.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 degrees) if dt=None, dt will be set such that the number of theta values is equal to the height 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
exception abel.tools.polar.CythonExtensionsNotBuilt[source]

Bases: exceptions.Exception

abel.tools.symmetry module

abel.tools.symmetry.get_image_quadrants(IM, reorient=True, symmetry_axis=None, use_quadrants=(True, True, True, True))[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 averaged. See Note.
  • use_quadrants (boolean tuple) – Include quadrant (Q0, Q1, Q2, Q3) in the symmetry combination(s) and final image
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 == averaged 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, odd_size=True, 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 |
    +--------+--------+
    
  • odd_size (boolean) – Whether final image is odd or even pixel size odd size will trim 1 row from Q1, Q0, and 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=False, dr=1, dt=None)[source]

Angular integration of the image.

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

Parameters:
  • IM (2D np.array) – The data image.
  • origin (tuple) – Image center coordinate relative to bottom-left corner defaults to rows//2+rows%2,cols//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).
  • dt (float) – Theta coordinate grid spacing in degrees. 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 np.array) – radial coordinates
  • speeds (1D np.array) – Integrated intensity array (vs radius).

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

Intensity variation in the angular coordinate, theta.

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

(optionally and more useful) returning intensity vs angle for defined radial ranges.

Parameters:
  • IM (2D np.array) – Image data
  • radial_ranges (list of tuples) – [(r0, r1), (r2, r3), ...] Evaluate the intensity vs angle for the radial ranges r0_r1, r2_r3, etc.
Returns:

  • intensity_vs_theta (2D np.array) – Intensity vs angle distribution for each selected radial range.
  • theta (1D np.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 = xs_total/4pi [ 1 + beta P2(cos theta) ]     Eq. (1)

where P2(x)=(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 np.array) – Angle coordinates, referenced to the vertical direction.
  • intensity (1D np.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
Returns:

  • (beta, error_beta) (tuple of floats) – The anisotropy parameters and the errors associated with each one.
  • (amplitude, error_amplitude) (tuple of floats) – Amplitude of signal and an error for each amplitude. Compare this with the data to check the fit.

abel.benchmark module

class abel.benchmark.AbelTiming(n=[301, 501], n_max_bs=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 global validity can be checked with array.all()
  • 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 abel.analytical, initialized) –
  • recon (1D ndarray) – a reconstruction (i.e. inverse abel) given by some PyAbel implementation
abel.benchmark.main()[source]

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.