pynx.utils: useful functions for fft, prime factor, pattern simulation, plotting

Array

pynx.utils.array.array_derivative(a, dr, mask=None, phase=False)

Compute the derivative of an array along a given direction

Parameters:
  • a – the complex array for which the derivative will be calculated

  • dr – the shift in pixels (with a value along each of the array dimensions) the value returned will be (a(r+dr)-a[r-dr])/2/norm2(dr), or if one of the values is masked, e.g. (a(r+dr)-a[r])/norm2(dr)

  • phase – if True, will return instead the derivative of the phase of the supplied complex array. Will return an error if the supplied array is not complex.

Returns:

a masked array of the gradient.

pynx.utils.array.center_array_2d(a, other_arrays=None, threshold=0.2, roi=None, iz=None)

Center an array in 2D so that its absolute value barycenter is in the middle. If the array is 3D, it is summed along the first axis to determine the barycenter, and all frames along the first axis are shifted. The array is ‘rolled’ so that values shifted from the right appear on the left, etc… Shifts are integer - no interpolation is done.

Parameters:
  • a – the array to be shifted, can be a floating point or complex 2D or 3D array.

  • other_arrays – can be another array or a list of arrays to be shifted by the same amount as a

  • threshold – only the pixels above the maximum amplitude * threshold will be used for the barycenter

  • roi – tuple of (x0, x1, y0, y1) corners coordinate of ROI to calculate center of mass

  • iz – if a.ndim==3, the centering will be done based on the center of mass of the absolute value summed over all 2D stacks. If iz is given, the center of mass will be calculated just on that stack

Returns:

the shifted array if only one is given or a tuple of the shifted arrays.

pynx.utils.array.crop(a: ndarray, margin=None, margin_f=None, shift=False)

Crop an array by removing values on each side

Parameters:
  • a – the array to crop

  • margin

    the margin to crop. This can either be:

    • an integer number (same number of pixels removed on each side and dimension)

    • a list/tuple of integers, one for each dimension. The same margin is applied on both sides

    • a list/tuple of integers, two for each dimension. A different margin can be applied on both sides.

  • margin_f – instead of giving the margin on each side as a number of pixels, it is possible to give it as a integer factor by which the final size will be divided, e.g. if margin_f=2, the size is divided by 2, and the margin is equal to 1/4 of the original size, for each dimension. This can be set either as an integer, with the same factor for each dimension, or as a list with one value for each dimension. Ignored if margin is given.

  • shift – if True, the input array is assumed to be fft-shifted (centred on the first array element), and will be shifted for cropping, and fft-shifted again for the output.

Returns:

the cropped array

pynx.utils.array.crop_around_support(obj: ndarray, sup: ndarray, margin=0)
Parameters:
  • obj – the array to be cropped (2D or 3D)

  • sup – the support, either a boolean or integer array, with the same dimensions as a, 0 (False) indicating the pixels outside the support.

  • margin – the number or pixels to be added on all sides of the array

Returns:

a tuple with (cropped array, cropped support), keeping only pixels inside the support

pynx.utils.array.fourier_shift(a, shift, axes=None, positivity=False)

Sub-pixel shift of an array. The return type will be the same as the input.

Parameters:
  • a – the array to shift, with N dimensions

  • shift – the shift along each axis

  • axes – a tuple of the axes along which the shift is performed. If None, all axes are transformed

  • positivity – if True, all values <0 will be set to zero for the output

Returns:

the fft-shifted array

pynx.utils.array.interp_linear(a: ndarray, dr, mask: ndarray | None = None, return_weight=False, use_numba=True)

Perform a bi/tri-linear interpolation of a 2D/3D array.

Parameters:
  • a – the 2D or 3D array to interpolate

  • dr – the shift (along each dimension) of the array for the interpolation

  • mask – the mask of values (True or >0 values are masked) to compute the interpolation

  • return_weight – if true, also return the weight of the interpolation, which should be equal to 1 if no point needed for the interpolation were masked, and a value between 0 and 1 otherwise.

Returns:

a masked interpolated array

pynx.utils.array.pad(a: ndarray, padding=None, padding_f=None, stack=False, value=0, shift=False)
Parameters:
  • a – the array to pad (2D, 2D stack, 3D)

  • padding – the number of pixels to add on each border. if this is an integer, the same margin is added on each border. This can be a list/tuple with a size the number of dimensions, so each dimension has a different margin on each side. Finally, there can be twice the number of dimensions so each side of each dimension uses a different margin.

  • padding_f – instead of giving a number of pixels to add, it is possible to give a factor by which the array size will be multiplied. This can either be an integer for all dimentions, or one value for eachi dimension. Ignored if padding is given.

  • stack – if True and the array is 3D, pad it as a stack of 2D arrays. In this case, padding should only contain the values for the two dimensions.

  • value – the value with which to fill the padded values

  • shift – if True, the input array is assumed to be fft-shifted, and will be centred before padding, and fft-shifted again upon return

Returns:

the padded array

pynx.utils.array.pad2(d, pad_width, mode='constant', mask2neg=False, **kwargs)

Pad a 2D array with an interface similar to numpy.pad(), with additional padding modes to avoid a discontinuity across the padded array border.

TODO: allow the same pad_width options as numpy.pad(), i.e. either a single value

or a list of tuples (before_i, after_i) for each dimension.

Parameters:
  • d – the 2D array to pad

  • pad_width – padding width, either a single value or one for each axis (pady, padx)

  • mode – the padding mode, either one from numpy, or among ‘reflect_linear’, ‘reflect_erf’, ‘reflect_sine’. The padded areas are equal to the reflected array, but with an interpolation allowing to avoid discontinuity at the array borders.

  • mask2neg – if True, the padded values are stored as -1-I_pad

  • kwargs – keywords passed to numpy.pad() in case it is used, ignored otherwise.

Returns:

the padded array.

pynx.utils.array.rebin(a, rebin_f, scale='sum', mask=None, mask_iobs_cdi=False, **kwargs)

Rebin a 2 or 3-dimensional array. If its dimensions are not a multiple of rebin_f, the array will be cropped.

Parameters:
  • a – the array to resize, which can also be a masked array

  • rebin_f – the rebin factor - pixels will be summed by groups of rebin_f x rebin_f (x rebin_f). This can also be a tuple/list of rebin values along each axis, e.g. rebin_f=(4,1,2) for a 3D array Instead of summing/averaging the pixels over the rebin box, it is also possible to select a sub-pixel by giving the shift for each dimension, e.g. with “rebin=4,1,2,0,0,1”, the extracted array will be a[0::4,0::1,1::2]

  • scale – if “sum” (the default), the array total will be kept. If “average”, the average pixel value will be kept. If “square”, the array is scaled so that (abs(a)**2).sum() is kept

  • mask – an array of values to be masked (0 or False, valid entries, >0 or True, should be masked). Alternatively, a can be a masked array.

  • mask_iobs_cdi – if True, special negative Iobs values used in pynx.cdi will be correctly handled.

Returns:

the array after rebinning. A masked array if mask is not None.

pynx.utils.array.upsample(a: ndarray, bin_f, scale='sum', interp=1)

Inverse operation of rebin, for a 2 or 3-dimensional array. The array dimensions are multiplied by rebin_f along each dimension. The added values are linearly interpolated, assuming circular periodicity across boundaries.

Parameters:
  • a – the array to resize

  • bin_f – the bin factor - each pixel/voxel will be transformed in bin_f x bin_f (x bin_f) pixels. This can also be a tuple/list of bin values along each axis, e.g. bin_f=(4,1,2) for a 3D array

  • scale – if “sum” (the default), the array total will be kept. If “average”, the average pixel value will be kept. If “square”, (abs(a)**2).sum() is kept

  • interp – if 1, a linear interpolation is used (which can be very slow). If 0, no interpolation.

Returns: the new array with dimensions multiplied by bin_f

Benchmark

pynx.utils.benchmark.benchmark_fft(device, do_plot=True, ndim=2, batch=1, max_size=4096, step_size=16, verbose=True)

Calculate FFT speed for one or several CUDA or OpenCL device(s)

Parameters:
  • device – the device to be tested. This can be a pycuda or pyopencl device, or a string to match against available devices names. It can also be a list/tuple of the two pycuda/pyopencl devices

  • do_plot – if True, plot using matplotlib

  • ndim – the number of dimensions for the FFT - either 2 or 3

  • batch – the batch number for the FFT,e.g. to perform 2D FFT on a 3D stack of ‘batch’ depth

  • max_size – the maximum (inclusive) size (along each FFT dimension) of the FFT. benchmark will stop increasing stop when an exception (likely memory allocation) is caught.

Returns:

a dictionary with the name of the device + CLFFT or CUFFT as key, and two vectors with fft size and gflops

Fourier shell correlation

class pynx.utils.fourier_shell_correlation.FSCPlot(img1, img2, snrt=[0.2071, 0.5], ring_thick=0, rad_apod=60, axial_apod=20, pixel_size=None)

Plot the FSC and threshold curves

Parameters:
  • img1 – first image (2D or 3D)

  • img2 – second image (2D or 3D)

  • snrt – power SNR for threshold computation. Options: SNRt = 0.5 -> 1 bit threshold for average SNRt = 0.2071 -> 1/2 bit threshold for average

  • ring_thick – thickness (in pixel units) of the frequency rings. Normally the pixels get assigned to the closest integer pixel ring in Fourier Domain. With ring_thick, each ring gets more pixels and more statistics.

  • rad_apod – radial apodisation width

  • axial_apod – axial apodisation width

plot(save_plot=None, plot_images=False, cmap=None, figsize=None)
Parameters:
  • save_plot – if given (string), save the plot to the given file name

  • plot_images – if True, the images will be plotted along the FSC plot

  • cmap – the colormap to be used. Useful for phase maps

  • figsize – thf igure size to use. If None, a default is used.

Returns:

Nothing. A dictionary with the plotted coordinates for all curves and coordinates is saved as self.d

class pynx.utils.fourier_shell_correlation.FourierShellCorr(img1, img2, snrt=0.2071, ring_thick=0, rad_apod=60, axial_apod=20)

Computes the Fourier Shell Correlation between image1 and image2, and computes the threshold funcion T of 1 or 1/2 bit. It can handle non-cube arrays, but it assumes that the voxel is isotropic. It applies a Hanning window of the size of the data to the data before the Fourier transform calculations to attenuate the border effects.

Reference: M. van Heel, M. Schatz, “Fourier shell correlation threshold criteria”, Journal of Structural Biology 151, 250-262 (2005) https://doi.org/10.1016/j.jsb.2005.05.009

@author: Julio Cesar da Silva (jdasilva@esrf.fr)

Parameters:
  • img1 – first image (2D or 3D)

  • img2 – second image (2D or 3D)

  • snrt – power SNR for threshold computation. Options: SNRt = 0.5 -> 1 bit threshold for average SNRt = 0.2071 -> 1/2 bit threshold for average

  • ring_thick – thickness (in pixel units) of the frequency rings. Normally the pixels get assigned to the closest integer pixel ring in Fourier Domain. With ring_thick, each ring gets more pixels and more statistics.

  • rad_apod – radial apodisation width

  • axial_apod – axial apodisation width

apodization()

Compute the Hanning window of the size of the data for the apodization

fouriercorr()

Compute FSC and threshold

nyquist()

Evaluate the Nyquist Frequency

ringthickness()

Define ring_thick

transverse_apodization()

Compute the Hanning window of the size of the data for the apodization

History

class pynx.utils.history.History

Class to record optimization history. It is used to store the parameters like the algorithm, negative log-likelihood (llk), chi^2, resolution, cputime, walltime as a function of the cycle number. Not all values need be stored for all cycles. The default fields initialized as python OrderedDict (the key being the cycle number) are:

  • ‘time’: records timeit.default_timer()

  • ‘epoch’: records time.time()

as_numpy_record_array(*args)

Get stored values for one or several keys (e.g. ‘time’, ‘llk’, ‘nb_photons’) as a numpy record array. The first entry is always the cycle number. If entries are missing for a given key and a cycle number, it is replaced by the next recorded value.

Parameters:

args – all the desired keys, e.g.: ‘time’, ‘epoch’, ‘llk’, ‘algorithm’… If no args are given, all available keys are returned.

Returns:

numpy record array, see https://docs.scipy.org/doc/numpy/user/basics.rec.html#record-arrays, or None if no history has been recorded

insert(cycle, **kwargs)

Store new values. if keys do not already exist, they are automatically added. ‘time’ and ‘epoch’ keys are automatically recorded, so need not be supplied. :param cycle: the current cycle :param **kwargs: e.g. llk=2e4, algorithm=’ML-Poisson’

Returns:

Nothing.

class pynx.utils.history.PynxOrderedDict

OrderedDict with easy access to the last value.

as_numpy_record_array(title='data')

Return dictionary as a numpy record array. Strings are converted to ASCII for h5py compatibility :return: the numpy record array, with two entries (cycle, value) per position

Math

pynx.utils.math.full_width(x, y, ratio=0.5, outer=False)

Determine the full-width from an XY dataset.

Parameters:
  • x – the abscissa (1d array) of the data to be fitted. Assumed to be in monotonic ascending order

  • y – the y-data as a 1d array. The absolute value of the array will be taken if it is complex.

  • ratio – the fraction (default=0.5 for FWHM) at which the width should be measured

  • outer – if True, the width will be measured by taking the outermost points which fall below the ratio to the maiximum. Otherwise, the width will be taken as the width around the maximum (regardless of secondary peak which may be above maw*ratio)

Returns:

the full width

pynx.utils.math.llk_euclidian(iobs, imodel)

Compute the Eucldian log-likelihood for a calculated intensity, given observed values. The value computed is normalized so that its asymptotic value (for a large number of points) is equal to the number of observed points. This model is valid if obs and calc are reasonably close.

Parameters:
  • iobs – the observed intensities

  • imodel – the calculated/model intensity

Returns:

The negative log-likelihood.

pynx.utils.math.llk_gaussian(iobs, imodel)

Compute the Gaussian log-likelihood for a calculated intensity, given observed values. The value computed is normalized so that its asymptotic value (for a large number of points) is equal to the number of observed points.

Parameters:
  • iobs – the observed intensities

  • imodel – the calculated/model intensity

Returns:

The negative log-likelihood.

pynx.utils.math.llk_poisson(iobs, imodel, imodel_min=0.1)

Compute the Poisson log-likelihood for a calculated intensity, given observed values. The value computed is normalized so that its asymptotic value (for a large number of points) is equal to the number of observed points

Parameters:
  • iobs – the observed intensities

  • imodel – the calculated/model intensity

  • imodel_min – the minimum accepted value for imodel (to avoid infinite llk when imodel=0 and iobs>0)

Returns:

The negative log-likelihood.

pynx.utils.math.ortho_modes(m, method='eig', verbose=False, return_matrix=False, nb_mode=None, return_weights=False)

Orthogonalize modes from a N+1 dimensional array or a list/tuple of N-dimensional arrays. The decomposition is such that the total intensity (i.e. (abs(m)**2).sum()) is conserved, and the modes are orthogonal, i.e. (mo[i]*mo[j].conj()).sum()=0 for i!=j

Parameters:
  • m – the stack of modes to orthogonalize along the first dimension.

  • method – either ‘eig’ to use eigenvalue decomposition, or ‘svd’ to use singular value decomposition.

  • verbose – if True, the matrix of coefficients will be printed

  • return_matrix – if True, return the linear combination matrix

  • nb_mode – the maximum number of modes to be returned. If None, all are returned.

  • return_weights – if True, will also return the relative weights of all the modes. This is useful if nb_mode is used, and only a partial list of modes is returned.

Returns:

an array (mo) with the same shape as given in input, but with orthogonal modes sorted by decreasing norm. Then if return_matrix is True, the matrix of linear coefficients is returned Then if return_weights is True, an array with the

pynx.utils.math.primes(n)

Returns the prime decomposition of n as a list

pynx.utils.math.smaller_primes(n, maxprime=13, required_dividers=(4,), decrease=True)

Find the closest integer <= or >=n (or list/array of integers), for which the largest prime divider is <=maxprime, and has to include some dividers. The default values for maxprime is the largest integer accepted by the clFFT library for OpenCL GPU FFT.

Parameters:
  • n – the integer number

  • maxprime – the largest prime factor acceptable

  • required_dividers – a list of required dividers for the returned integer.

  • decrease – if True (thed default), the integer returned will be <=n, otherwise it will be >=n.

Returns:

the integer (or list/array of integers) fulfilling the requirements

pynx.utils.math.test_smaller_primes(n, maxprime=13, required_dividers=(4,))

Test if the largest prime divider is <=maxprime, and optionally includes some dividers.

Parameters:
  • n – the integer number for which the prime decomposition will be tested

  • maxprime – the maximum acceptable prime number. This defaults to the largest integer accepted by the clFFT library for OpenCL GPU FFT.

  • required_dividers – list of required dividers in the prime decomposition. If None, this test is skipped.

Returns:

True if the conditions are met.

Pattern simulation

pynx.utils.pattern.fibonacci_urchin(dsize=256, nb_rays=100, r_max=64, nb_rings=8)

Calculate a binary urchin (in 3D).

Parameters:
  • dsize – size in pixels for the 2D array with the star data

  • nb_rays – number of radial branches. Must be > 0

  • r_max – maximum radius in pixels

  • nb_rings – number of rings (the rays will have some holes between successive rings)

Returns:

a 3D array with the binary urchin.

pynx.utils.pattern.get_img(index=0)

Returns image (numpy array) from scipy.misc

Parameters:

index – 0-> return scipy.misc.face()[:,:,1], cropped to 512x512 1 -> return 512x512 scipy.misc.ascent()

pynx.utils.pattern.siemens_star(dsize=512, nb_rays=36, r_max=None, nb_rings=8, cheese_holes_nb=0, cheese_hole_max_radius=5, cheese_hole_spiral_period=0)

Calculate a binary Siemens star.

Parameters:
  • dsize – size in pixels for the 2D array with the star data

  • nb_rays – number of radial branches for the star. Must be > 0

  • r_max – maximum radius for the star in pixels. If None, dsize/2 is used

  • nb_rings – number of rings (the rays will have some holes between successive rings)

  • cheese_holes_nb – number of cheese holes other the entire area, resulting more varied frequencies

  • cheese_hole_max_radius – maximum axial radius for the holes (with random radius along x and y). If the value is negative, the radius is fixed instead of random.

  • cheese_hole_spiral_period – instead of randomly distributing holes, giving an integer N number for this parameter will generate holes located on an Archimedes spiral, every N pixels along the spiral, which also has a step of N pixels. The pattern of holes is aperiodic.

Returns:

a 2D array with the Siemens star.

pynx.utils.pattern.spiral_archimedes(a, n)

“ Creates np points spiral of step a, with a between successive points on the spiral. Returns the x,y coordinates of the spiral points.

This is an Archimedes spiral. the equation is:

r=(a/2*pi)*theta the stepsize (radial distance between successive passes) is a the curved absciss is: s(theta)=(a/2*pi)*integral[t=0->theta](sqrt(1*t**2))dt

pynx.utils.pattern.spiral_fermat(dmax, n)

” Creates a Fermat spiral with n points distributed in a circular area with diameter<= dmax. Returns the x,y coordinates of the spiral points. The average distance between points can be roughly estimated as 0.5*dmax/(sqrt(n/pi))

http://en.wikipedia.org/wiki/Fermat%27s_spiral

Phase

pynx.utils.phase.minimize_grad_phase(im, mask_thr=0.3, center_phase=None, global_min=False, mask=None, rebin_f=None)

Minimises the gradient in the phase of the input image im.

Parameters:
  • im – the input complex 2D image for which the phase gradient needs to be minimized

  • mask_thr – only pixel amplitudes above mask_thr*abs(im).max() will be taken into account

  • center_phase – after minimizing the gradient, center the phase around this value (in radians). pi/3 gives a nice contrast

  • global_min – if True after a quick local minimization of the gradient, the phase differences will be minimized over the whole object. Use this for flat-phased objects

  • mask – a 2d mask array (True = masked points) giving the pixels which should be excluded from the optimization. This is cumulative with mask_thr

  • rebin_f – an integer>1 can be given to compute the phase gradient on a rebinned array for faster evaluation.

Returns:

a tuple with (phase corrected image, correction array, mask, linear coefficients)

pynx.utils.phase.phase_tilt(im, grad=(0, 0), offset=0)

Multiplies the phase of the complex image by the gradient and adds phase offset. DEPRECATED. scipy.ndimage.fourier_shift should be used instead

Parameters:
  • im – input image

  • grad – gradient of the phase as a multiple of 2*pi across the image

  • offset – offset of the phase as a multiple of 2*pi

pynx.utils.phase.remove_phase_ramp(d: ndarray, mask_thr=None, center_phase=None, mask=None, niter=2)

Removes a linear ramp in the phase of the input array (2D or 3D). This uses scikit-image phase_cross_correlation in Fourier space for a fast evaluation, by calculating the shift of FT(d) compared to FT(abs(d)). 3D correction requires scikit-image>=0.15.0.

WARNING: This works efficiently on a pure phase ramp, with some noise, but cannot be used in general - e.g. on an object with a phase vortex, FT(d) and FT(abs(d)) can be quite different and cannot easily be correlated.

Parameters:
  • d – the data array for which the phase ramp must be removed.

  • mask_thr – only pixel amplitudes above mask_thr*abs(im).max() will be taken into account

  • center_phase – after minimizing the gradient, center the phase around this value (in radians). pi/3 gives a nice contrast.

  • mask – a mask array (True = masked points) giving the pixels which should be excluded from the optimization. This is cumulative with mask_thr. Masked pixels are corrected, but not used to evaluate the phase gradient.

  • niter – number of iterations to perform to remove ramp d (default=2)

Returns:

a tuple with (phase corrected array, linear coefficients). The linear coefficients correspond to the total amplitude of the ramp along the given axis, in radians.

pynx.utils.phase.shift_phase_zero(obj, percent=5, mask=None, origin=0, stack=True, verbose=False)

Shift the phase of the given object, so that the phase range begins at 0 radians, and hopefully avoid wrapped phases. If the percentile range from (percent, 100-percent) cannot be lowered to less than 5.5 radians minus twice percent*2*pi/199, it is assumed that the object has phase wrapping and no correction can be made, so the object is returned unchanged.

Parameters:
  • obj – the complex object for which a shift of the phase will be calculated. If the object is 3D, the phase shift is evaluated only on the first mode obj[0].

  • percent – the range of the phase will be evaluated using np.percentile from ‘percent’ to ‘100-percent’. This is used to avoid the influence of noise

  • mask – if given, only the non-masked area (for which mask > 0) will be taken into account

  • origin – the desired origin of the phase. 0 by default, but e.g. -3 can be used for a [-pi;pi] display

  • stack – if True, treat the object as a stack of 2D object (e.g. the object modes in ptychography), and calculate the phase shift on the first layer, but apply to all layers. if True, the mask should be 2D

Returns:

the object corrected by a constant phase shift, i.e. obj * exp(1j * dphi)

pynx.utils.phase.unwrap_phase(d: array, method='skimage', medfilt_n=3)

Extract an unwrap the phase from a complex data array. Currently implemented for 2D only.

Parameters:
  • d – the complex data array to be unwrapped

  • method – either ‘skimage’ to use skimage.restoration.unwrap_phase, or ‘gradient-x’ to use the normalised complex array gradient to determine the smooth phase gradient and then integrate along the x-axis (last dimension), or ‘gradient-y’ to use a gradient and final integration along y.

  • medfilt_n – size of the median kernel filtering on the initial array or phase. This is only used to evaluate the phase wrapping - the final returned phase array is not filtered.

Returns:

the unwrapped phase in radians.

Phase retrieval transfer function

pynx.utils.phase_retrieval_transfer_function.nyquist(shape)

Evaluate the Nyquist Frequency

pynx.utils.phase_retrieval_transfer_function.plot_prtf(freq, fnyquist, prtf, iobs_shell=None, nbiobs_shell=None, pixel_size=None, title=None, file_name=None, fig_num=101)

Plot the phase retrieval transfer function

Parameters:
  • freq – the frequencies for which the phase retrieval transfer function was calculated

  • fnyquist – the nyquist frequency

  • prtf – the phase retrieval transfer function

  • iobs_shell – the integrated observed intensity per frequency shell (same size as prtf), which will be plotted against a second y-axis and can give some information about the amount of scattering.

  • nbiobs_shell – number of points per iobs_shell (to be able to compute the average Iobs per pixel)

  • pixel_size – the pixel size in metres, for resolution axis

  • title – the figure title

  • file_name – if given, the plot will be saved to this file (’.png’ or ‘.pdf’)

pynx.utils.phase_retrieval_transfer_function.prtf(icalc, iobs, mask=None, ring_thick=5, shell_averaging_method='before', norm_percentile=100, scale=False)

Compute the phase retrieval transfer function, given calculated and observed intensity. Note that this function assumes that calc and obs are orthonormal arrays.

Parameters:
  • icalc – the calculated intensity array (origin at array center), either 2D or 3D

  • iobs – the observed intensity (origin at array center)

  • mask – the mask, values > 0 indicating bad pixels in the observed intensity array

  • ring_thick – the thickness of each shell or ring, in pixel units

  • shell_averaging_method – by default (‘before’) the amplitudes are averaged over the shell before the PRTF=<calc>/<obs> is computed. If ‘after’ is given, then the ratio of calc and observed amplitudes will first be calculated (excluding zero-valued observed pixels), and then averaged to compute the PRTF=<calc/obs>.

  • norm_percentile – the output PRTF is normalised so that the ‘norm_percentile’ percentile value is 1. If None, no normalisation is performed.

  • scale – if True (default=False), the calculated array is scaled so that the non-masked integrated intensity are equal. Should only be used with norm_percentile=None

Returns:

a tuple with the (frequency, frequency_nyquist, prtf, iobs), where frequency, prtf are masked arrays, fnyquist the Nyquist frequency, and iobs the integrated observed intensity per shell

pynx.utils.phase_retrieval_transfer_function.ring_thickness(shape)

Define ring_thick

Plot utils

pynx.utils.plot_utils.colorwheel(text_col='black', fs=16, ax=None)

Color wheel for phases in hsv colormap.

Parameters:
  • text_col – colour of text

  • fs – fontsize in points

  • fig – an instance of matplotlib.figure.Figure (for offline rendering) or None

Returns:

Nothing. Displays a colorwheel in the current or supplied figure.

pynx.utils.plot_utils.complex2rgbalin(s, gamma=1.0, smax=None, smin=None, percentile=(None, None), alpha=(0, 1), type='uint8')

Returns RGB image with with colour-coded phase and linear amplitude in brightness. Optional exponent gamma is applied to the amplitude.

Parameters:
  • s – the complex data array (likely 2D, but can have higher dimensions)

  • gamma – gamma parameter to change the brightness curve

  • smax – maximum value (brightness = 1). If not supplied and percentile is not set, the maximum amplitude of the array is used.

  • smin – minimum value(brightness = 0). If not supplied and percentile is not set, the maximum amplitude of the array is used.

  • percentile – a tuple of two values (percent_min, percent_max) setting the percentile (between 0 and 100): the smax and smin values will be set as the percentile value in the array (see numpy.percentile). These two values (when not None) supersede smax and smin. Example: percentile=(0,99) to scale the brightness to 0-1 between the 1% and 99% percentile of the data amplitude.

  • alpha – the minimum and maximum value for the alpha channel, normally (0,1). Useful to have different max/min alpha when going through slices of one object

  • type – either ‘float’: values are in the [0..1] range, or ‘uint8’ (0..255) (new default)

Returns:

the RGBA array, with the same diemensions as the input array, plus one additional R/G/B/A dimension appended.

pynx.utils.plot_utils.complex2rgbalin_dark(s, gamma=1, smin=None, smax=None, percentile=None, phase_shift=0)

Returns RGB image based on a HSV projection with hue=phase, saturation=1 and value =amplitude. This yields a black-background image contrary to complex2rgbalin

Parameters:
  • s – the 2D complex data array

  • gamma – gamma parameter to change the brightness curve

  • smax – maximum value (value = 1). If not supplied and percentile is not set, the maximum amplitude is used.

  • smin – minimum value(value = 0). If not supplied and percentile is not set, the minimum amplitude is used.

  • percentile – a tuple of two values (percent_min, percent_max) setting the percentile (between 0 and 100): the smax and smin values will be set as the percentile value in the array (see numpy.percentile). These two values (when not None) supersede smax and smin. Example: percentile=(0,99) to scale the brightness to 0-1 between the 1% and 99% percentile of the array amplitude.

  • phase_shift – shift the phase origin by this amount (in radians)

Returns:

the RGB array, with values between 0 and 1.

pynx.utils.plot_utils.complex2rgbalog(s, amin=0.5, dlogs=2, smax=None, type='uint8')

Returns a 2D RGBA image with colour-coded phases and log10(amplitude) in brightness.

Parameters:
  • s – the 2D complex data array to be converted to RGB

  • amin – the minimum value for the alpha channel

  • dlogs – amplitude range displayed, in log10 units

  • smax – if specified, all values above max will be displayed with an alpha value of 1

  • type – either ‘float’: values are in the [0..1] range, or ‘uint8’ (0..255) (new default)

Returns:

pynx.utils.plot_utils.insertColorwheel(left=0.7, bottom=0.15, width=0.1, height=0.1, text_col='black', fs=16)

Inserts color wheel to the current axis.

pynx.utils.plot_utils.phase2rgb(s)

Crates RGB image with colour-coded phase from a complex array.

Parameters:

s – a complex numpy array

Returns:

an RGBA numpy array, with one additional dimension added

pynx.utils.plot_utils.showCplx(im, mask=0, pixSize_um=1, showGrid=True, amplitudeLog=False, maskPhase=False, maskPhaseThr=0.01, cmapAmplitude='gray', cmapPhase=<matplotlib.colors.LinearSegmentedColormap object>, scalePhaseImg=True, suptit=None, fontSize=20, suptit_fontSize=10, show_what='amplitude_phase', hideTicks=False)

Displays AMPLITUDE_PHASE or REAL_IMAG (‘show_what’) of the complex image in two subfigures.

pynx.utils.plot_utils.show_obj_probe(obj, probe, tit1='Object', tit2='Probe', stit=None, fig_num=100, pixel_size_object=None, scan_area_obj=None, scan_area_probe=None, scan_pos=None, invert_yaxis=False, **kwargs)

Live plot of the object and probe phase and amplitude. Only the first modes are shown.

Parameters:
  • obj – the object to plot

  • probe – the probe to plot

  • tit1 – title for the object

  • tit2 – title for the probe

  • stit – suptitle for the figure (e.g. the current algorithm, llk,…)

  • fig_num – the figure number to use, so plotting always appears in the same one

  • pixel_size_object – pixel size for the object

  • scan_area_obj – the 2D array of the object scan area (1 inside, 0 outside), which will be sued to mask the object

  • scan_pos – a tuple (x,y) of the scan positions in meters.

  • invert_yaxis – if True, the plot will be inverted along the Y axis- this is used for near field ptychography so the view corresponds to the full field view, rather than using the motor axes as reference.