Common
- ptychography40.reconstruction.common.apply_matrix(sources, matrix, target_shape)[source]
Apply a transformation matrix generated by
image_transformation_matrix()
to a stack of images.- Parameters:
sources (array-like) – Array of shape (n, sy, sx) where (sy, sx) is the
source_shape
parameter ofimage_transformation_matrix()
.matrix (array-like) – Matrix generated by
image_transformation_matrix()
or equivalenttarget_shape (Tuple[int, int]) –
source_shape
parameter ofimage_transformation_matrix()
. The result will be reshaped to(n, ) + target_shape
- ptychography40.reconstruction.common.diffraction_to_detector(lamb, diffraction_shape, pixel_size_real, pixel_size_detector, cy, cx, flip_y=False, scan_rotation=0.0)[source]
Generate a function that transforms pixel coordinates from diffraction of real space to detector coordinates.
When performing a forward calculation where a wave front is passed through an object function and projected in the far field, the projection is a Fourier transform. Each pixel in the fourier-transformed slice of real space corresponds to a diffraction angle in radian. This angle is a function of real space dimensions and wavelength.
For ptychography, this forward-projected beam is put in relation with detector data. The projection and detector, however, are usually not calibrated in such a way that each pixel corresponds exactly to one diffracted pixel/beam from the object and illumination function. This function applies the correct scale, rotation and handedness to match the forward-projected beam with the detector data, given the necessary input parameters.
cy, cx, flip_y and scan_rotation are chosen to correspond to the parameters in :meth:libertem.api.Context.create_com_analysis`.
This function is designed to create a
affine_transformation
parameter forimage_transformation_matrix()
.- Parameters:
lamb (float) – Wavelength in m
diffraction_shape (Tuple[int, int]) – Shape of the diffracted area
pixel_size_real (float or tuple or ndarray) – Pixel size in m for the diffracted shape. Can be a Tuple[float, float] or array with two values for y and x
pixel_size_detector (float or tuple or ndarray) – Pixel size in radian of the detector. Can be a Tuple[float, float] or array with two values for y and x. For free propagation into the far-field, this is the detector pixel size in m divided by the camera length for small angles.
cy (float) – Y position of the central beam on the detector in pixel.
cx (float) – X position of the central beam on the detector in pixel.
flip_y (bool) – Flip the y axis of the detector coordinates
scan_rotation (float) – Scan rotation in degrees
- Returns:
transform – A function that accepts pixel coordinates in diffracted space as an array of shape (n, 2) with (y, x). Upper left corner is (0, 0). It returns pixel coordinates on the detector as floats of shape (n, 2).
- Return type:
callable(coords : numpy.ndarray) -> numpy.ndarray
- ptychography40.reconstruction.common.fftshift_coords(reconstruct_shape)[source]
Generate a function that performs an FFT shift of coordinates.
On the detector, the central beam is near the center, while for native FFT the center is at the (0, 0) position of the result array. Instead of fft-shifting the result of a projection calculation to match the detector directly, one can instead calculate a transformation that picks from an unshifted FFT result in such a way that it matches a diffraction pattern in its native layout. This allows to combine the FFT shift with additional transformation steps. See also
image_transformation_matrix()
.- Parameters:
- Returns:
Callable(coords – Function that accepts target coordinates with shape (n, 2), kind int and calculates the source coordinates with shape (n, 2), kind int so that taking from the source coordinates and placing into the target coordinates performs an FFT shift.
- Return type:
- ptychography40.reconstruction.common.get_shifted(arr_shape: tuple, tile_origin: tuple, tile_shape: tuple, shift: tuple)[source]
Calculate the slices to cut out a shifted part of a 2D source array and place it into a target array, including tiling support.
This works with negative and positive integer shifts.
- ptychography40.reconstruction.common.ifftshift_coords(reconstruct_shape)[source]
Generate a function that performs an inverse FFT shift of coordinates.
On the detector, the central beam is near the center, while for native FFT the center is at the (0, 0) position of the result array. Instead of fft-shifting the result of a projection calculation to match the detector, one can instead calculate a transformation that picks from the detector data in such a way that an inverse FFT shift is performed. This allows to combine the inverse FFT shift with additional transformations. See also
image_transformation_matrix()
.- Parameters:
- Returns:
Callable(coords – Function that accepts target coordinates with shape (n, 2), kind int and calculates the source coordinates with shape (n, 2), kind int so that taking from the source coordinates and placing into the target coordinates performs an inverse FFT shift.
- Return type:
- ptychography40.reconstruction.common.image_transformation_matrix(source_shape, target_shape, affine_transformation, pre_transform=None, post_transform=None)[source]
Construct a sparse matrix that transforms a flattened source image stack to a flattened target image stack.
A sparse matrix prodct can be a highly efficient method to apply a set of transformations to an image in one pass. This function constructs a sparse matrix that picks values from a source image to fill each pixel of the target image by first applying
pre_transform()
to the target image indices to map them into a 2D vector space, then projecting the pixel outline in this vector space usingaffine_transformation()
to calculate the source pixel coordinates, and then usingpost_transform()
to map the source coordinates to indices in the source image. If the projected pixel is of size 1.5 or smaller in the source coordinates, the closest integer is chosen. If it is larger, the average of pixels within the projected pixel outline is chosen. This corresponds to scaling with order=0 inscipy.ndimage.zoom()
.pre_transform() and :code:`post_transform()
can also be used for shifting the center ofaffine_transformation()
.- Parameters:
source_shape (Tuple[int, int]) – Shape of source and target image for bounds checking and index raveling
target_shape (Tuple[int, int]) – Shape of source and target image for bounds checking and index raveling
affine_transformation (callable(coords -> coords)) – Transformation that maps intermediate coordinates, i.e. the result of applying
pre_transform()
, to float source coordinates. It should be continuous, strictly monotone and approximately affine for the size of one pixel.diffraction_to_detector()
can be used to generate a suitable coordinate transformation function.pre_transform (callable(coords) -> coords) – Map target image indices to coordinates, typically euclidean.
pre_transform()
should not change the scale of the coordinates. It is designed to be something likeifftshift_coords()
to un-scramble target coordinates so that coordinates that are close in the source image are also close in the un-scrambled intermediate coordinates generated by this function. This is identity by default.post_transform (callable(coords) -> coords(int)) – Map source image coordinates, typically euclidean, to source image indices.
post_transform()
should not change the scale of the coordinates. By default it isnp.round(...).astype(int)
.
- Returns:
Shape np.prod(source_shape), np.prod(target_shape)
- Return type:
scipy.sparse.csc_matrix
- ptychography40.reconstruction.common.rolled_object_aggregation_cpu(obj_out, updates, shifts, fftshift=False)[source]
Aggregate shifted updates to an object function
This function accumulates updates that are shifted relative to the object function using addition and rolls the indices within the object if necesssary. Optionally, it can fftshift the updates while integrating. Doing this in one loop allows to reduce the number of calls for each shift and reduces overall memory transfer.
- Parameters:
obj_out (numpy.ndarray) – 2D array with the object, modified in-place by this function
updates (numpy.ndarray) – Array with updates, shape (n, …)
shifts (numpy.ndarray) – Array with shift vectors, shape (n, 2), kind int
fftshift (bool) – Read the updates fft-shifted from
updates
- ptychography40.reconstruction.common.rolled_object_aggregation_cuda(obj_out, updates, shifts, fftshift=False)[source]
Numba CUDA version of
rolled_object_aggregation_cpu()
- ptychography40.reconstruction.common.rolled_object_probe_product_cpu(obj, probe, shifts, result_out, ifftshift=False)[source]
Multiply object and shifted illumination
This function combines several steps that are relevant for ptychographic reconstruction:
Multiply an object function with a shifted illumination
Roll the object function indices around the edges
Optionally, perform an ifftshift to prepare the data for subsequent FFT
These steps are combined in a single loop since each requires significant memory transfer if they are performed step-by-step. For performance reasons it doesn’t perform a free subpixel shift, but picks the best option from a set of pre-calculated shifted versions.
See
shifted_probes()
for a function to calculate the shifted versions.- Parameters:
obj (numpy.ndarray) – 2D array with the object
probe (numpy.ndarray) – 4D array with subpixel shifts of the probe, last two dimensions same size or smaller than obj.
shifts (numpy.ndarray) – Array with shift vectors, shape (n, 2), kind float
result_out (numpy.ndarray) – Array where the result is placed. Shape (n, ) + probe.shape
ifftshift (bool) – place the product ifft-shifted into
result_out
- Returns:
subpixel_indices – The first two indices for
probe
- Return type:
np.ndarray
- ptychography40.reconstruction.common.rolled_object_probe_product_cuda(obj, probe, shifts, result_out, ifftshift=False)[source]
Numba CUDA version of
rolled_object_probe_product_cpu()
- ptychography40.reconstruction.common.shifted_probes(probe, bins)[source]
Calculated subpixel-shifted versions of the probe
- Parameters:
probe (numpy.ndarray) –
bins (int or Tuple[int, int]) – Number of antialiasing steps in y and x axis. Can be int as well
- Returns:
probes – 4D, shape bins + probe.shape or (bins, bins) + probe.shape if bins is an int
- Return type: