operators#
Defines implementations for all operators.
All of the solvers, rely on operators including forward and adjoint operators. In tike, forward and adjoint operators are paired as fwd and adj methods of an Operator.
In this way, multiple solvers (e.g. ePIE, gradient descent, SIRT) implemented in Python can share the same core operators and can be upgraded to better operators in the future.
All operator methods accept the array type that matches the output of their asarray() method.
A base class for Operators. |
|
An alignment operator composed of pad, flow, and rotate operations. |
|
A Laminography operator. |
|
Provides a multi-plan per-device cache for CuPy FFT. |
|
A 2D Convolution operator with linear interpolation. |
|
Map input 2D arrays to new coordinates by Lanczos interpolation. |
|
A Laminography operator. |
|
Pad a stack of 2D images to the same shape but with unique pad_widths. |
|
Extract (zero-padded) patches from images at provided positions. |
|
A Fourier-based free-space propagation using CuPy. |
|
A Ptychography operator. |
|
Rotate a stack of 2D images along last two dimensions. |
|
Shift last two dimensions of an array using Fourier method. |
- class tike.operators.Alignment[source]#
An alignment operator composed of pad, flow, and rotate operations.
The operations are applied in the aforementioned order.
Please see the help for the Pad, Flow, and Rotate operations for description of arguments.
- __new__(**kwargs)#
- class tike.operators.Bucket[source]#
A Laminography operator.
Laminography operators to simulate propagation of the beam through the object for a defined tilt angle. An object rotates around its own vertical axis, nz, and the beam illuminates the object some tilt angle off this axis.
- n#
The pixel width of the cubic reconstructed grid.
- Type
int
- tilt#
The tilt angle; the angle between the rotation axis of the object and the light source. π / 2 for conventional tomography. 0 for a beam path along the rotation axis.
- Type
float32
- Parameters
u ((nz, n, n) complex64) – The complex refractive index of the object. nz is the axis corresponding to the rotation axis.
data ((ntheta, n, n) complex64) – The complex projection data of the object.
theta (array-like float32) – The projection angles; rotation around the vertical axis of the object.
- __new__(**kwargs)#
- class tike.operators.CachedFFT[source]#
Provides a multi-plan per-device cache for CuPy FFT.
A class which inherits from this class gains the _fft2, _fftn, and _ifft2 methods which provide automatic plan caching for the CuPy FFTs.
This plan cache differs from the cache included in CuPy>=8 because it is NOT per-thread. This allows us to use threadpool.map() and allows us to destroy the cache manually.
- __init__()#
- __new__(**kwargs)#
- class tike.operators.Convolution[source]#
A 2D Convolution operator with linear interpolation.
Compute the product two arrays at specific relative positions.
- nscan#
The number of scan positions at each angular view.
- Type
int
- probe_shape#
The pixel width and height of the (square) probe illumination.
- Type
int
- nz, n
The pixel width and height of the reconstructed grid.
- Type
int
- ntheta#
The number of angular partitions of the data.
- Type
int
- Parameters
psi ((..., nz, n) complex64) – The complex wavefront modulation of the object.
probe (complex64) – The (…, nscan, nprobe, probe_shape, probe_shape) or (…, 1, nprobe, probe_shape, probe_shape) complex illumination function.
nearplane (complex64) – The (…., nscan, nprobe, probe_shape, probe_shape) wavefronts after exiting the object.
scan ((..., nscan, 2) float32) – Coordinates of the minimum corner of the probe grid for each measurement in the coordinate system of psi. Vertical coordinates first, horizontal coordinates second.
- __new__(**kwargs)#
- class tike.operators.Flow[source]#
Map input 2D arrays to new coordinates by Lanczos interpolation.
Uses Lanczos interpolation for a non-affine deformation of a series of 2D images.
- __init__()#
- __new__(**kwargs)#
- class tike.operators.Lamino[source]#
A Laminography operator.
Laminography operators to simulate propagation of the beam through the object for a defined tilt angle. An object rotates around its own vertical axis, nz, and the beam illuminates the object some tilt angle off this axis.
- n#
The pixel width of the cubic reconstructed grid.
- Type
int
- tilt#
The tilt angle; the angle between the rotation axis of the object and the light source. π / 2 for conventional tomography. 0 for a beam path along the rotation axis.
- Type
float32 [radians]
- Parameters
u ((nz, n, n) complex64) – The complex refractive index of the object. nz is the axis corresponding to the rotation axis.
data ((ntheta, n, n) complex64) – The complex projection data of the object.
theta (array-like float32 [radians]) – The projection angles; rotation around the vertical axis of the object.
- __new__(**kwargs)#
- class tike.operators.Operator[source]#
A base class for Operators.
An Operator is a context manager which provides the basic functions (forward and adjoint) required solve an inverse problem.
Operators may be composed into other operators and inherited from to provide additional implementations to the ones provided in this library.
- __init__()#
- __new__(**kwargs)#
- class tike.operators.Pad[source]#
Pad a stack of 2D images to the same shape but with unique pad_widths.
By default, no padding is applied and/or the padding is applied symmetrically.
- __init__()#
- __new__(**kwargs)#
- class tike.operators.Patch[source]#
Extract (zero-padded) patches from images at provided positions.
- Parameters
images ((..., H, W) complex64) – The complex wavefront modulation of the object.
positions ((..., N, 2) float32) – Coordinates of the minimum corner of the patches in the image grid.
patches ((..., N * nrepeat, width+, width+) complex64) – OR (…, L, width+, width+) complex64 The extracted (zero-padded) patches. (N * nrepeat) = K * L, K >= nrepeat
patch_width (int) – The width of the unpadded patches.
- __init__()#
- __new__(**kwargs)#
- class tike.operators.Propagation[source]#
A Fourier-based free-space propagation using CuPy.
Take an (…, N, N) array and apply the Fourier transform to the last two dimensions.
- detector_shape#
The pixel width and height of the nearplane and farplane waves.
- Type
int
- cost#
The function to be minimized when solving a problem.
- Type
(data-like, farplane-like) -> float
- grad#
The gradient of cost.
- Type
(data-like, farplane-like) -> farplane-like
- Parameters
nearplane ((..., detector_shape, detector_shape) complex64) – The wavefronts after exiting the object.
farplane ((..., detector_shape, detector_shape) complex64) – The wavefronts hitting the detector respectively. Shape for cost functions and gradients is (nscan, 1, 1, detector_shape, detector_shape).
Changed in version 0.25.0: Removed the model parameter and the cost(), grad() functions. Use the cost and gradient functions directly instead.
- __init__(detector_shape, norm='ortho', **kwargs)[source]#
- Parameters
detector_shape (int) –
norm (str) –
- __new__(**kwargs)#
- class tike.operators.Ptycho[source]#
A Ptychography operator.
Compose a diffraction and propagation operator to simulate the interaction of an illumination wavefront with an object followed by the propagation of the wavefront to a detector plane.
- Parameters
detector_shape (int) – The pixel width and height of the (square) detector grid.
nz (int) – The pixel width and height of the reconstructed grid.
n (int) – The pixel width and height of the reconstructed grid.
probe_shape (int) – The pixel width and height of the (square) probe illumination.
propagation (
Operator
) – The wave propagation operator being used.diffraction (
Operator
) – The object probe interaction operator being used.data ((..., FRAME, WIDE, HIGH) float32) – The intensity (square of the absolute value) of the propagated wavefront; i.e. what the detector records.
farplane ((..., POSI, 1, SHARED, detector_shape, detector_shape) complex64) – The wavefronts hitting the detector respectively.
probe ({(..., 1, 1, SHARED, WIDE, HIGH), (..., POSI, 1, SHARED, WIDE, HIGH)} complex64) – The complex illumination function.
psi ((..., WIDE, HIGH) complex64) – The wavefront modulation coefficients of the object.
scan ((..., POSI, 2) float32) – Coordinates of the minimum corner of the probe grid for each measurement in the coordinate system of psi. Coordinate order consistent with WIDE, HIGH order.
Changed in version 0.25.0: Removed the model and ntheta parameters.
- __init__(detector_shape, probe_shape, nz, n, propagation=<class 'tike.operators.cupy.propagation.Propagation'>, diffraction=<class 'tike.operators.cupy.convolution.Convolution'>, norm='ortho', **kwargs)[source]#
Please see help(Ptycho) for more info.
- Parameters
detector_shape (int) –
probe_shape (int) –
nz (int) –
n (int) –
propagation (Type[Propagation]) –
diffraction (Type[Convolution]) –
norm (str) –
- __new__(**kwargs)#
- class tike.operators.Rotate[source]#
Rotate a stack of 2D images along last two dimensions.
- Parameters
angle (float) – The desired rotation in radians. Operation skipped if angle is None.
cval (complex64) – The value to use for filling regions that rotated from outside the original image.
- __init__()#
- __new__(**kwargs)#
- class tike.operators.Shift[source]#
Shift last two dimensions of an array using Fourier method.
- __init__()#
- __new__(**kwargs)#
- tike.operators.gaussian(data, intensity)[source]#
The Gaussian model objective function.
- Parameters
data ((N, M, M)) – The measured diffraction data
intensity ((N, M, M)) – The modeled intensity
- Return type
float
- tike.operators.gaussian_each_pattern(data, intensity)[source]#
The Gaussian model objective function per diffraction pattern.
- Parameters
data ((N, M, M)) – The measured diffraction data
intensity ((N, M, M)) – The modeled intensity
- Returns
costs – The objective function for each pattern.
- Return type
(N, )
- tike.operators.gaussian_grad(data, farplane, intensity)[source]#
The gradient of the Gaussian model objective function
- Parameters
data ((N, M, M)) – The measured diffraction data
intensity ((N, M, M)) – The modeled intensity
farplane ((N, K, L, M, M)) –
- Return type
cupy.ndarray
- tike.operators.poisson(data, intensity)[source]#
The Poisson model objective function.
- Parameters
data ((N, M, M)) – The measured diffraction data
intensity ((N, M, M)) – The modeled intensity
- Return type
float
- tike.operators.poisson_each_pattern(data, intensity)[source]#
The Poisson model objective function per diffraction pattern.
- Parameters
data ((N, M, M)) – The measured diffraction data
intensity ((N, M, M)) – The modeled intensity
- Returns
costs – The objective function for each pattern.
- Return type
(N, )