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.

Module for operators utilizing the CuPy library.

This module implements the forward and adjoint operators using CuPy. This removes the need for interface layers like pybind11 or SWIG because kernel launches and memory management may by accessed from Python.

Operator

A base class for Operators.

Alignment

An alignment operator composed of pad, flow, and rotate operations.

CachedFFT

Provides a multi-plan cache for CuPy FFT.

Convolution

A 2D Convolution operator with linear interpolation.

Flow

Map input 2D arrays to new coordinates by Lanczos interpolation.

Lamino

A Laminography operator.

Pad

Pad a stack of 2D images to the same shape but with unique pad_widths.

Propagation

A Fourier-based free-space propagation using CuPy.

Ptycho

A Ptychography operator.

Rotate

Rotate a stack of 2D images along last two dimensions.

Shift

Shift last two dimensions of an array using Fourier method.

class tike.operators.cupy.Alignment[source]

Bases: tike.operators.cupy.operator.Operator

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.

adj(rotated, flow, shift, unpadded_shape, angle, padded_shape=None, cval=0.0)[source]

Perform the adjoint operator.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
fwd(unpadded, shift, flow, padded_shape, angle, unpadded_shape=None, cval=0.0)[source]

Perform the forward operator.

inv(rotated, flow, shift, unpadded_shape, angle, padded_shape=None, cval=0.0)[source]
class tike.operators.cupy.Convolution(probe_shape, nz, n, ntheta=None, detector_shape=None, **kwargs)[source]

Bases: tike.operators.cupy.operator.Operator

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.

adj(nearplane, scan, probe, psi=None, overwrite=False)[source]

Combine probe shaped patches into a psi shaped grid by addition.

adj_probe(nearplane, scan, psi, overwrite=False)[source]

Combine probe shaped patches into a probe.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
fwd(psi, scan, probe)[source]

Extract probe shaped patches from the psi at each scan position.

The patches within the bounds of psi are linearly interpolated, and indices outside the bounds of psi are not allowed.

class tike.operators.cupy.Flow[source]

Bases: tike.operators.cupy.operator.Operator

Map input 2D arrays to new coordinates by Lanczos interpolation.

Uses Lanczos interpolation for a non-affine deformation of a series of 2D images.

adj(g, flow, filter_size=5, cval=0.0)[source]

Remap individual pixels of f with Lanczos filtering.

Parameters
  • g ((…, H, W) complex64) – A stack of deformed arrays.

  • flow ((…, H, W, 2) float32) – The displacements to be applied to each pixel along the last two dimensions. Operation skipped when flow is None.

  • filter_size (int) – The width of the Lanczos filter. Automatically rounded up to an odd positive integer.

  • cval (complex64) – This value is used for interpolation from points outside the grid.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
fwd(f, flow, filter_size=5, cval=0.0)[source]

Remap individual pixels of f with Lanczos filtering.

Parameters
  • f ((…, H, W) complex64) – A stack of arrays to be deformed.

  • flow ((…, H, W, 2) float32) – The displacements to be applied to each pixel along the last two dimensions. Operation skipped when flow is None.

  • filter_size (int) – The width of the Lanczos filter. Automatically rounded up to an odd positive integer.

  • cval (complex64) – This value is used for interpolation from points outside the grid.

inv(g, flow, filter_size=5, cval=0.0)[source]
class tike.operators.cupy.Lamino(n, tilt, eps=0.001, **kwargs)[source]

Bases: tike.operators.cupy.cache.CachedFFT, tike.operators.cupy.operator.Operator

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.

adj(data, theta, overwrite=False, **kwargs)[source]

Perform the adjoint Laminography transform.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
cost(data, theta, obj)[source]

Cost function for the least-squres laminography problem

fwd(u, theta, **kwargs)[source]

Perform the forward Laminography transform.

gather(Fe, x, n, m, mu)[source]
grad(data, theta, obj)[source]

Gradient for the least-squares laminography problem

scatter(f, x, n, m, mu)[source]
class tike.operators.cupy.Operator[source]

Bases: abc.ABC

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.

adj(**kwargs)[source]

Perform the adjoint operator.

classmethod asarray(*args, device=None, **kwargs)[source]
classmethod asnumpy(*args, **kwargs)[source]
fwd(**kwargs)[source]

Perform the forward operator.

class tike.operators.cupy.Pad[source]

Bases: tike.operators.cupy.operator.Operator

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.

adj(padded, corner=None, unpadded_shape=None, **kwargs)[source]

Strip the edges from the padded images.

Parameters
  • corner ((N, 2)) – The min corner of the images in the padded array.

  • padded_shape (3-tuple) – The desired shape after padding. First element should be N.

  • unpadded_shape (3-tuple) – See padded_shape.

  • cval (complex64) – The value to use for padding.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
fwd(unpadded, corner=None, padded_shape=None, cval=0.0, **kwargs)[source]

Pad the unpadded images with cval.

Parameters
  • corner ((N, 2)) – The min corner of the images in the padded array.

  • padded_shape (3-tuple) – The desired shape after padding. First element should be N.

  • unpadded_shape (3-tuple) – See padded_shape.

  • cval (complex64) – The value to use for padding.

inv(padded, corner=None, unpadded_shape=None, **kwargs)

Strip the edges from the padded images.

Parameters
  • corner ((N, 2)) – The min corner of the images in the padded array.

  • padded_shape (3-tuple) – The desired shape after padding. First element should be N.

  • unpadded_shape (3-tuple) – See padded_shape.

  • cval (complex64) – The value to use for padding.

class tike.operators.cupy.Patch[source]

Bases: tike.operators.cupy.operator.Operator

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) – The extracted (zero-padded) patches.

  • patch_width (int) – The width of the unpadded patches.

adj(positions, patches, images=None, patch_width=None, height=None, width=None, nrepeat=1)[source]

Perform the adjoint operator.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
fwd(images, positions, patches=None, patch_width=None, height=None, width=None, nrepeat=1)[source]

Perform the forward operator.

class tike.operators.cupy.Propagation(detector_shape, model='gaussian', **kwargs)[source]

Bases: tike.operators.cupy.cache.CachedFFT, tike.operators.cupy.operator.Operator

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

model

The type of noise model to use for the cost functions.

Type

string

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 (ntheta, nscan, 1, 1, detector_shape, detector_shape).

  • data, intensity ((ntheta, nscan, detector_shape, detector_shape) complex64) – data is the square of the absolute value of farplane. data is the intensity of the farplane.

adj(farplane, overwrite=False, **kwargs)[source]

Adjoint Fourier-based free-space propagation operator.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
fwd(nearplane, overwrite=False, **kwargs)[source]

Forward Fourier-based free-space propagation operator.

class tike.operators.cupy.Ptycho(detector_shape, probe_shape, nz, n, ntheta=1, model='gaussian', propagation=<class 'tike.operators.cupy.propagation.Propagation'>, diffraction=<class 'tike.operators.cupy.convolution.Convolution'>, **kwargs)[source]

Bases: tike.operators.cupy.operator.Operator

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

  • model (string) – The type of noise model to use for the cost functions.

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

adj(farplane, probe, scan, psi=None, overwrite=False, **kwargs)[source]

Please see help(Ptycho) for more info.

adj_probe(farplane, scan, psi, overwrite=False, **kwargs)[source]

Please see help(Ptycho) for more info.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
cost(data, psi, scan, probe)float[source]

Please see help(Ptycho) for more info.

fwd(probe, scan, psi, **kwargs)[source]

Please see help(Ptycho) for more info.

grad_probe(data, psi, scan, probe, mode=None)[source]

Compute the gradient with respect to the probe(s).

Parameters

mode (list(int)) – Only return the gradient with resepect to these probes.

grad_psi(data, psi, scan, probe)[source]

Please see help(Ptycho) for more info.

class tike.operators.cupy.Rotate[source]

Bases: tike.operators.cupy.operator.Operator

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.

adj(rotated, angle, cval=0.0)[source]

Perform the adjoint operator.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
fwd(unrotated, angle, cval=0.0)[source]

Perform the forward operator.

inv(rotated, angle, cval=0.0)[source]
class tike.operators.cupy.Shift[source]

Bases: tike.operators.cupy.cache.CachedFFT, tike.operators.cupy.operator.Operator

Shift last two dimensions of an array using Fourier method.

adj(a, shift, overwrite=False, cval=None)[source]

Perform the adjoint operator.

classmethod asarray(*args, device=None, **kwargs)
classmethod asnumpy(*args, **kwargs)
fwd(a, shift, overwrite=False, cval=None)[source]

Apply shifts along last two dimensions of a.

Parameters
  • array (…, H, W) float32 – The array to be shifted.

  • shift (…, 2) float32 – The the shifts to be applied along the last two axes.

inv(a, shift, overwrite=False, cval=None)

Perform the adjoint operator.