reborn.misc package¶
Submodules¶
reborn.misc.correlate module¶
- reborn.misc.correlate.plan_masked_pearson_cc_2d(x_mask, y_mask, max_shift: (<class 'float'>, None) = None, min_unmasked: int = 100) dict [source]¶
Create a plan to speed up the masked_pearson_cc_2d function. Determines the minimal slice of the image that is needed for the correlation. Stores the sliced masks and FFTs of the sliced/padded masks to avoid re-calculation.
- Parameters:
x (
ndarray
) – First imagex_mask (
ndarray
) – First image masky (
ndarray
) – Second imagey_mask (
ndarray
) – Second image maskmax_shift (float) – Set the correlations to zero where the shift (distance, in pixel units) is larger than this number. If the value is less than one, it is interpreted as a fraction of the radius of the smallest edge of the image.
min_unmasked (int) – Set the correlations to zero where the number of jointly unmasked pixels is less than this
- Returns:
dict
- reborn.misc.correlate.masked_pearson_cc_2d(x, x_mask, y, y_mask, max_shift=None, min_unmasked=9, plan=None)[source]¶
Cross correlate two masked arrays. Unlike an ordinary cross correlation, this cross correlation normalizes the result so that it spans from -1 to 1 as defined by the Pearson correlation formula. Masked pixels are handled by excluding them from the cross correlations; i.e. for every shift in the cross correlation, only the pixels for which both pixels are unmasked will contribute to the cross correlation. FFT operations are used for speed. Options are provided to set correlations to zero wherever the number of jointly unmasked pixels is too low (poor statistics) or where the shift is too large (avoids FFT artifacts).
The returned Pearson cross correlation is shifted such that the shifts of an NxM array are as follows:
nv, nw = image.shape vshifts = np.arange(nv*2) - nv wshifts = np.arange(nw*2) - nw
- Parameters:
x (
ndarray
) – First imagex_mask (
ndarray
) – First image masky (
ndarray
) – Second imagey_mask (
ndarray
) – Second image maskmax_shift (float) – Set the correlations to zero where the shift (distance, in pixel units) is larger than this number. If the value is less than one, it is interpreted as a fraction of the radius of the smallest edge of the image.
min_unmasked (int) – Set the correlations to zero where the number of jointly unmasked pixels is less than this
- Returns:
(
ndarray
) – The Pearson correlation
reborn.misc.interpolate module¶
- reborn.misc.interpolate.trilinear_interpolation(densities, vectors, corners=None, deltas=None, x_min=None, x_max=None, out=None, strict_types=True)[source]¶
Perform a trilinear interpolation on a 3D array. An arbitrary set of sample points in the form of an \(N\times 3\) array may be specified.
Notes
This function behaves as if the density is periodic; points that lie out of bounds will wrap around. This might change in the future, in which case a keyword argument will be added so that you may explicitly decide what to do in the case of points that lie outside of the grid. Note that periodic boundaries avoid the need for conditional statements within a for loop, which probably makes the function faster. For now, if you think you have points that lie outside of the grid, consider handling them separately.
You may specify the output array, which is useful if you wish to simply add to an existing array that you have already allocated. This can make your code faster and reduce memory. Beware: the out array is not over-written – the underlying fortran function will add to the existing out array.
Only double precision arrays (both real and complex are allowed) at the fortran level. You may pass in other types, but they will be converted to double (or double complex) before passing to the fortran function.
Make sure that all your arrays are c-contiguous.
An older version of this code allowed the arguments corners and deltas. They are discouraged because we aim to standardize on the x_min and x_max arguments documented below. They may be removed in the future.
The shape of the 3D array is inferred from the densities argument.
- Parameters:
densities (numpy array) – A 3D density array.
vectors (numpy array) – An Nx3 array of vectors that specify the points to be interpolated.
x_min (float or numpy array) – A 3-element vector specifying the center of the corner voxel of the 3D array. If a float is passed instead, it will be replicated to make a 3D array.
x_max (float or numpy array) – Same as x_min, but specifies the opposite corner, with larger values than x_min.
out (numpy array) – If you don’t want the output array to be created (e.g for speed), provide it here.
strict_types (bool) – Set this to False if you don’t mind your code being slow due to the need to convert datatypes (i.e. copy arrays) on every function call. Default: True.
- Returns:
numpy array
- reborn.misc.interpolate.trilinear_insertion(densities, weights, vectors, insert_vals, corners=None, deltas=None, x_min=None, x_max=None)[source]¶
Perform the “inverse” of a trilinear interpolation . That is, take an arbitrary set of sample values along with their 3D vector locations and “insert” them into a 3D grid of densities. The values are distributed amongst the nearest 8 grid points so that they sum to the original insert value.
FIXME: Be more clear about the mathematical operation that this function performs…
Notes
This function behaves as if the density is periodic; points that lie out of bounds will wrap around. This might change in the future, in which case a keyword argument will be added so that you may explicitly decide what to do in the case of points that lie outside of the grid. Note that periodic boundaries avoid the need for conditional statements within a for loop, which probably makes the function faster. For now, if you think you have points that lie outside of the grid, consider handling them separately.
You may specify the output array, which is useful if you wish to simply add to an existing 3D array that you have already allocated. This can make your code faster and reduce memory. Beware: the out array is not over-written – the underlying fortran function will add to the existing
densities
array.Only double precision arrays (both real and complex are allowed).
Make sure that all your arrays are c-contiguous.
An older version of this code allowed the arguments
corners
anddeltas
. They are discouraged because we aim to standardize on thex_min
andx_max
arguments documented below. They may be removed in the future.The shape of the 3D array is inferred from the
densities
argument.
- Parameters:
densities (numpy array) – A 3D array containing the densities, into which values are inserted. Note that an “insertion” means that the
insert_vals
below are multiplied byweights
below.weights (numpy array) – A 3D array containing weights. These are needed in order to perform a weighted average. After calling this function one or more times, the average densities are calculated by dividing
densities
byweights
. Be mindful of divide-by-zero errors.vectors (numpy array) – The 3D vector positions corresponding to the values to be inserted.
insert_vals (numpy array) – The values to be inserted into the 3D map. They are multiplied by weights before being inserted into the densities map.
x_min (float or numpy array) – A 3-element vector specifying the center of the corner voxel of the 3D array. If a float is passed instead, it will be replicated to make a 3D array.
x_max (float or numpy array) – Same as x_min, but specifies the opposite corner, with larger values than x_min.
- Returns:
None. This function modifies the densities and weights arrays; it returns nothing.
- reborn.misc.interpolate.trilinear_insertion_factor(densities, weight_factor, vectors, insert_vals, corners=None, deltas=None, x_min=None, x_max=None)[source]¶
Performs trilinear insert with a factor being multiplied onto the weights.
Notes
This function behaves as if the density is periodic; points that lie out of bounds will wrap around. This might change in the future, in which case a keyword argument will be added so that you may explicitly decide what to do in the case of points that lie outside of the grid. Note that periodic boundaries avoid the need for conditional statements within a for loop, which probably makes the function faster. For now, if you think you have points that lie outside of the grid, consider handling them separately.
You may specify the output array, which is useful if you wish to simply add to an existing 3D array that you have already allocated. This can make your code faster and reduce memory. Beware: the out array is not over-written – the underlying fortran function will add to the existing
densities[1,:,:,:]
array.Make sure that all your arrays are c-contiguous.
An older version of this code allowed the arguments
corners
anddeltas
. They are discouraged because we aim to standardize on thex_min
andx_max
arguments documented below. They may be removed in the future.The shape of the 3D array is inferred from the
densities
argument.
- Parameters:
densities (numpy array) – A 4D array containing the densities, into which values are inserted. Note that an “insertion” means that the
insert_vals
below are multiplied byweights
below.weight_factor (numpy array) – A number that gets multiplied with the trilinear insertion weights.
vectors (numpy array) – The 3D vector positions corresponding to the values to be inserted.
insert_vals (numpy array) – The values to be inserted into the 3D map. They are multiplied by weights before being inserted into the densities map.
x_min (float or numpy array) – A 3-element vector specifying the center of the corner voxel of the 3D array. If a float is passed instead, it will be replicated to make a 3D array.
x_max (float or numpy array) – Same as x_min, but specifies the opposite corner, with larger values than x_min.
- Returns:
None. This function modifies the densities and weights arrays; it returns nothing.
- reborn.misc.interpolate.trilinear_insertions(densities, vectors, insert_vals, corners=None, deltas=None, x_min=None, x_max=None)[source]¶
Performs multiple trilinear inserts
Notes
This function behaves as if the density is periodic; points that lie out of bounds will wrap around. This might change in the future, in which case a keyword argument will be added so that you may explicitly decide what to do in the case of points that lie outside of the grid. Note that periodic boundaries avoid the need for conditional statements within a for loop, which probably makes the function faster. For now, if you think you have points that lie outside of the grid, consider handling them separately.
You may specify the output array, which is useful if you wish to simply add to an existing 3D array that you have already allocated. This can make your code faster and reduce memory. Beware: the out array is not over-written – the underlying fortran function will add to the existing
densities[1,:,:,:]
array.Make sure that all your arrays are c-contiguous.
An older version of this code allowed the arguments
corners
anddeltas
. They are discouraged because we aim to standardize on thex_min
andx_max
arguments documented below. They may be removed in the future.The shape of the 3D array is inferred from the
densities
argument.
- Parameters:
densities (numpy array) – A 4D array containing the densities, into which values are inserted.
vectors (numpy array) – The 3D vector positions corresponding to the values to be inserted.
insert_vals (numpy array) – The values to be inserted into the 3D map. They are multiplied by weights before being inserted into the densities map.
x_min (float or numpy array) – A 3-element vector specifying the center of the corner voxel of the 3D array. If a float is passed instead, it will be replicated to make a 3D array.
x_max (float or numpy array) – Same as x_min, but specifies the opposite corner, with larger values than x_min.
- Returns:
None. This function modifies the densities array; it returns nothing.
- reborn.misc.interpolate.trilinear_insert(data_coord, data_val, x_min, x_max, n_bin, mask, boundary_mode='truncate', verbose=True)[source]¶
Trilinear insertion on a regular grid with arbitrarily positioned sample points.
This function returns two arrays, dataout and weightout. weightout is a 3D array containing the accumulated trilinear weights. dataout is the accumulated trilinearly inserted values. One needs to divide dataout by weightout (taking care to deal with zeros in weightout) to get the correct trilinear insertion result. This is so that the function can be used to sum over many trilinearly inserted arrays in for example a 3D diffracted intensity merge.
Note 1: All input arrays should be C contiguous. Note 2: This code will break if you put a 1 in any of the N_bin entries. Note 3: The boundary is defined as [x_min-0.5, x_max+0.5).
- Parameters:
data_coord (Nx3
ndarray
) – Coordinates (x,y,z) of the data points that you wish to insert into the regular grid.data_val (Nx1
ndarray
) – The values of the data points to be inserted into the grid.x_min (1x3
ndarray
) – (x_min, y_min, z_min)x_max (1x3
ndarray
) – (x_max, y_max, z_max)n_bin (1x3
ndarray
) – Number of bins in each direction (N_x, N_y, N_z)mask (Nx1
ndarray
) – Specify which data points to ignore. Non-zero means use, zero means ignore.boundary_mode (str) – Specify how the boundary should be treated. Options are: (1) “truncate” - Ignores all points outside the insertion volume. (2) “periodic” - Equivalent to wrapping around.
- Returns:
2-element tuple containing the following
dataout (3D |ndarray|) : Trilinearly summed values that needs to be divided by weightout to give the trilinearly inserted values.
weightout (3D |ndarray|) : Cumulative trilinear weights.
reborn.misc.polar module¶
- reborn.misc.polar.q_bin_indices(n_q_bins, q_bin_size, q_min, qs)[source]¶
Polar q bin indices. Given an arbitrary array of q magnitudes along with the binning specification (the center of the first bin, the bin size and the number of bins), output an array of bin indices. The indices correspond to Python indices, with the first index equal to zero.
- reborn.misc.polar.p_bin_indices(n_p_bins, p_bin_size, p_min, ps)[source]¶
Polar azimuthal (“phi”) bin indices. Given an arbitrary array of phis along with the binning specification (the center of the first bin, the bin size and the number of bins), output an array of bin indices. The indices correspond to Python indices, with the first index equal to zero.
- reborn.misc.polar.stats(data, q, p, weights=None, n_q_bins=100, q_min=0, q_max=30000000000.0, n_p_bins=360, p_min=0, p_max=3.141592653589793, sum_=None, sum2=None, w_sum=None)[source]¶
Gather polar-binned stats. Namely, the weighted sum, weighted sum of squares, and the sum of weights.
- Parameters:
data (
ndarray
)q
p
weights
n_q_bins
q_min
q_max
n_p_bins
p_min
p_max
sum
sum2
w_sum
Returns:
- class reborn.misc.polar.PolarPADAssembler(pad_geometry=None, beam=None, n_q_bins=50, q_range=None, n_phi_bins=None, phi_range=None)[source]¶
Bases:
object
A class for converting PAD data to polar coordinates.
- Parameters:
pad_geometry (
PADGeometryList
) – PAD Geometry.beam (
Beam
) – Beam information.n_q_bins (int) – Number of q bins.
q_range (tuple) – Minimum and maximum q bin centers. If None, the range is [0, maximum q in PAD].
n_phi_bins (int) – Number of phi bins.
phi_range (tuple) – Minimum and maximum phi bin centers. If None, the full 2*pi ring is assumed.
- classmethod from_resolution(pad_geometry, beam, sample_diameter, oversample=2)[source]¶
Instantiate PolarPADAssembler the correct way, with resolution and oversample ratio.
- Parameters:
pad_geometry (
PADGeometryList
) – PAD Geometry.beam (
Beam
) – Beam information.sample_diameter (float) – Size of sample (L).
oversample (int) – Oversmaple ratio (s=2 is default).
reborn.misc.rotate module¶
- reborn.misc.rotate.kabsch(A, A0)[source]¶
Finds the rotation matrix that will bring a set of vectors A into alignment with another set of vectors A0. Uses the Kabsch algorithm implemented in scipy.spatial.transform.Rotation.align_vectors
- class reborn.misc.rotate.Rotate3D(f3d, keep_last_even_k=False)[source]¶
Bases:
object
Base class to rotate a 3D array of double precision complex numbers in 3-dimensions. The function works by rotating each 2D sections of the 3D array via three shears, as described by Unser et al. (1995) “Convolution-based interpolation for fast, high-quality rotation of images.” IEEE Transactions on Image Processing, 4:1371.
- Note 1: The input array must be 3d, double complex, and have all three
dimension sizes equal. Otherwise, it raises a ValueError exception.
- Note 2: If you don’t want wrap-arounds, make sure the input array, f,
is zero-padded to at least sqrt(2) times the largest dimension of the desired object.
- Parameters:
f (3D |ndarray|) – The 3D input array. f is the corresponding
output. (class member for)
keep_last_even_k (bool) – default False. The last k for even N has an ambiguous
only (sign. Keeping one sign)
complex (makes real data become)
so
data. (this is always False for real)
- property f¶
- class reborn.misc.rotate.Rotate3DLinearShear(f3d, keep_last_even_k=False)[source]¶
Bases:
object
Same as Rotate3D but with the shears done with linear interpolation instead of FFTs. No _rotate3Dxr and _rotate3Dyr
- Note 1: The input array must be 3d, double complex, and have all three
dimension sizes equal. Otherwise, it raises a ValueError exception.
- Note 2: If you don’t want wrap-arounds, make sure the input array, f,
is zero-padded to at least sqrt(2) times the largest dimension of the desired object.
- Parameters:
f (3D |ndarray|) – The 3D input array. f is the corresponding
output. (class member for)
keep_last_even_k (bool) – default False. The last k for even N has an ambiguous
only (sign. Keeping one sign)
complex (makes real data become)
so
data. (this is always False for real)
- property f¶
- class reborn.misc.rotate.Rotate3Dlegacy(f3d)[source]¶
Bases:
Rotate3D
This is identical to rotate3D except that the shear orders are the same as in Joe’s original code, and the last lastk value is not zeroed. It is somewhat less efficient since the arrays are transposed and then transposed back.
- class reborn.misc.rotate.Rotate3Dvkfft(f3d, keep_last_even_k=False)[source]¶
Bases:
Rotate3D
This should give results identical to rotate3D for sizes that are products of powers of 2,3,5,7,11,13. It uses the pyvkfft wrapper to VkFFT to perform Fourier transforms on a gpu. Since this is my first opencl code, it is no doubt written inefficiently.
- property f¶
reborn.misc.spherical module¶
- class reborn.misc.spherical.ylmIntegration(L)[source]¶
Bases:
object
Class for calculating coefficients \(a_{lm}\) in the spherical harmonics expansion
(1)¶\[f(\theta, \phi) = \sum_{lm} a_{lm} Y_{lm}(\theta, \phi) ;.\]The spherical harmonics \(Y_{lm}(\theta, \phi)\) are the same as the scipy function
scipy.special.sph_harm
, HOWEVER, the \(\theta\) and \(\phi\) variables in the scipy documentation are switched when compared to the common physics convention used here. The normalization is defined so that setting \(f(\theta, \phi) = Y_{l'm'}(\theta, \phi)\) you will result in(2)¶\[a_{lm} = \delta(l-l')\delta(m-m')\]The initialization routine calculates and stores the angular spacing \(d\theta = \pi/(2L+2)\), the \(2L+2\) equally spaced midpoint \(\theta\) values (we exclude the poles \(\theta=0\) and \(\theta= \pi\)), the \(2L+2\) array of equally spaced \(\phi\) values (starting at \(\phi=0\)). Additional quantities (e.g. integration weights) are also calculated and stored.
- Parameters:
L (int) – Maximum \(L\) for \(Y^*_{lm} Y_{lm}\) to be correctly integrated.
- calc_ylmcoef(f)[source]¶
Calculates expansion coefficients \(a_{lm}\). The function \(f(\theta, \phi)\) must be provided as an
ndarray
of real or complex values corresponding to the stored grid of \(\theta\), \(\phi\) coordinates. The coefficients \(a_{lm}\) are returned in anndarray
of shape (\(L+1\), \(2L+1\)), with \(a_{lm} = 0\) for unused \(m\) values \(-L < m < L\). Numpy wrap-around indexing is used so that negative \(m\) array indices correspond to coefficients with negative indices.