Python Reference

speckle_tracking.make_mask(data, thresh=400, mask_gaps=False)[source]

Make a binary True/False mask from input data

Parameters
  • data (ndarray) – Input data, of shape (N, M, L).

  • thresh (float, optional) – The threshold value of the statistical measure below which pixels are considered to be good (mask == False).

Returns

mask – Boolean array of shape (M, L), where True indicates a good pixel and False a bad pixel.

Return type

ndarray

speckle_tracking.make_whitefield(data, mask, verbose=True)[source]

Estimate the image one would obtain without the sample in the beam.

This is done by taking the median value at each pixel along the first axis, then we try to fill in bad / zero pixels (where mask == False).

Parameters
  • data (ndarray) – Input data, of shape (N, M, L).

  • mask (ndarray) – Boolean array of shape (M, L), where True indicates a good pixel and False a bad pixel.

  • verbose (bool, optional) – print what I’m doing.

Returns

W – Float array of shape (M, L) containing the estimated whitefield.

Return type

ndarray

speckle_tracking.guess_roi(W, verbose=True)[source]

Find the rectangular region that contains most of the whitefield.

Parameters
  • W (ndarray) – The whitefield, that is the image one obtains without a sample in place.

  • verbose (bool, optional) – print what I’m doing.

Returns

roi – e.g. roi = [10, 400, 23, 500], indicates that most of the interesting data in a frame will be in the region: frame[roi[0]:roi[1], roi[2]:roi[3]]

Return type

list

speckle_tracking.fit_defocus(data, x_pixel_size, y_pixel_size, z, wav, mask, W, roi, verbose=True, **kwargs)[source]

Estimate the focus to sample distance.

This routine uses speckle_tracking.fit_thon_rings to estimate the defocus.

Parameters
  • data (ndarray) – Input data, of shape (N, M, L). Where N = number of frames M = number of pixels along the slow scan axis of the detector L = number of pixels along the fast scan axis of the detector

  • x_pixel_size (float) – The side length of a detector pixel in metres, along the slow scan axis.

  • y_pixel_size (float) – The side length of a detector pixel in metres, along the fast scan axis.

  • z (float) – The distance between the focus and the detector in meters.

  • wav (float) – The wavelength of the imaging radiation in metres.

  • W (ndarray) – The whitefield array of shape (M, L). This is the image one obtains without a sample in place.

  • roi (array_like) – Length 4 list of integers e.g. roi = [10, 400, 23, 500], indicates that most of the interesting data in a frame will be in the region: frame[roi[0]:roi[1], roi[2]:roi[3]]

  • verbose (bool, optional) – print what I’m doing.

  • kwargs (dict) – keyword arguments that are passed on to any functions called.

Returns

  • defocus (float) – The average focus to sample distance

  • dz (float) – The difference between the average defocus and the defocus along the fast and slow axes of the detector: z_ss = defocus + dz z_fs = defocus - dz defocus > 0 –> sample downstream of focus

  • res (dict) – Contains diagnostic information

See also

fit_thon_rings(), fit_defocus_registration()

speckle_tracking.fit_thon_rings(data, x_pixel_size, y_pixel_size, z, wav, mask, W, roi, centre=None, sig=10, edge_pix=5, window=30, rad_range=None, verbose=True, **kwargs)[source]

Find the focus to sample distance by fitting Thon rings to power spectrum.

This is done by generating a filtered power spectrum of the data. Then fitting concentric rings to this profile. The fitting parameters are then used to estimate the horizontal and vertical focus to sample distance.

Parameters
  • data (ndarray) – Input data, of shape (N, M, L). Where N = number of frames M = number of pixels along the slow scan axis of the detector L = number of pixels along the fast scan axis of the detector

  • x_pixel_size (float) – The side length of a detector pixel in metres, along the slow scan axis.

  • y_pixel_size (float) – The side length of a detector pixel in metres, along the fast scan axis.

  • z (float) – The distance between the focus and the detector in meters.

  • wav (float) – The wavelength of the imaging radiation in metres.

  • W (ndarray) – The whitefield array of shape (M, L). This is the image one obtains without a sample in place.

  • roi (array_like) – Length 4 list of integers e.g. roi = [10, 400, 23, 500], indicates that most of the interesting data in a frame will be in the region: frame[roi[0]:roi[1], roi[2]:roi[3]]

  • centre (array_like, optional) – Length 2 list of integers designating the centre of the Gaussian that is applied to the data before averaging. Default value is None in which case the centre of roi region is used.

  • edge_pix (integer, optional) – Number of edge pixels to mask in the power spectrum.

  • window (int, optional) – The sidelength of a square used to flatten the power spectrum contrast through the use of min / max filters.

  • rad_range (array_like or None, optional) – Length 2 list [min val., max val.] of the pixel radius within the power spectrum to fit the Thon rings. If None, then this is set to: [min(10, W.shape[0], W.shape[1]), max(W.shape[0], W.shape[1])//2]

  • verbose (bool, optional) – print what I’m doing.

Returns

  • defocus (float) – The average focus to sample distance

  • dz (float) – The difference between the average defocus and the defocus along the fast and slow axes of the detector: z_ss = defocus + dz z_fs = defocus - dz defocus > 0 –> sample downstream of focus

  • res (dict) – Contains diagnostic information:

    res[‘thon_display’]array_like, float

    shows the thon rings and the fit rings in one quadrant of the array.

    res[‘bd’]float

    The ratio of the real to the imaginary part of the refractive index, \(\beta_\lambda / \delta_\lambda\)

Notes

This routine fits the following function to the modulus of the power spectrum:

\[f(q_{ss}, q_{fs}) = p(q) | \delta_\lambda \sin(\pi\lambda z_2 q'^2) + \beta_\lambda \cos(\pi\lambda z_2 q'^2) |\]

where \(p(q)\) is a q dependent profile that depends on the details of the object, \(\lambda\) is the wavelength and q’ is given by:

\[q'^2 = (1+\frac{z_2}{z_{ss}})q_{ss}^2 + (1+\frac{z_2}{z_{fs}})q_{fs}^2\]

subject to:

\[z = \frac{1}{2}(z_{ss} + z_{fs}) + z_2\]

where \(z_{ss}\) and \(z_{fs}\) are the distance between the focal plane and the sample along the slow and fast scan axes respectively and \(z_2\) is the distance between the sample and the detector.

speckle_tracking.generate_pixel_map(shape, translations, basis, x_pixel_size, y_pixel_size, z, defocus_fs, defocus_ss=None, dss=None, dfs=None, verbose=True)[source]

Generate the pixel mapping based on the imaging geometry.

speckle_tracking.make_reference(data, mask, W, dij_n, pixel_map, roi=None, subpixel=False, verbose=True, minimum_overlap=None, sig=None)
Parameters
  • data (ndarray) – Input data, of shape (N, M, L).

  • mask (ndarray) – Boolean array of shape (M, L), where True indicates a good pixel and False a bad pixel.

  • W (ndarray) – Float array of shape (M, L) containing the estimated whitefield.

  • dij_n (ndarray) – Float array of shape (N, 2) containing the object translations that have been mapped onto the detector’s frame of reference.

  • pixel_map (ndarray, (2, M, L)) –

    An array containing the pixel mapping between a detector frame and the object, such that:

    \[\begin{split}I^{z_1}_{\phi}[n, i, j] = W[i, j] I^\infty[&\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0,\\ &\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0]\end{split}\]

  • subpixel (bool, optional) – If True then use bilinear subpixel interpolation non-integer pixel mappings.

  • verbose (bool, optional) – print what I’m doing.

  • minimum_overlap (float or None, optional) – Default is None. If float then the the object will be set to -1 where the number of data points contributing to that value is less than “minimum_overlap”.

Returns

  • I (ndarray) – Float array of shape (U, V), this is essentially an object map given by:

    \[\begin{split}I_0[i, j] &= \frac{\sum_n M[u_n, v_n] W[u_n, v_n] I^{z_1}_\phi[n, u_n, v_n]}{\sum_n M[u_n, v_n] W[u_n, v_n]^2 } \\\end{split}\]

    where:

    \[\begin{split}\begin{align} u_n[i,j] &= \text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0 \\ u_n[i,j] &= \text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0 \\ \end{align}\end{split}\]

    see Notes for more.

  • n0 (float) – Slow scan offset to the pixel mapping such that:

    \[\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0 \ge -0.5 \quad\text{for all } i,j\]
  • m0 (float) – Fast scan offset to the pixel mapping such that:

    \[\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0 \ge -0.5 \quad\text{for all } i,j\]

    -0.5 is chosen rather than 0 because integer coordinates are defined at the centre of the physical pixel locations.

Notes

\[I_{\phi, n}(\mathbf{x}) = W(\mathbf{x})I_0(\mathbf{x} - \frac{\lambda z}{2\pi} \nabla \Phi(x)-\Delta x_n, \bar{z}_\Phi)\]

\(M, W\) are the mask and whitefield arrays respectively. U and V are the pixel dimensions of \(I_0\) given by:

\[\begin{split}U &= \text{max}(\text{ij}_\text{map}[0, i, j]) - \text{min}(\Delta ij_n[0]) + n_0 \\ V &= \text{max}(\text{ij}_\text{map}[1, i, j]) - \text{min}(\Delta ij_n[1]) + m_0\end{split}\]
speckle_tracking.update_pixel_map(data, mask, W, O, pixel_map, n0, m0, dij_n, search_window=None, grid=None, roi=None, subpixel=False, subsample=1.0, interpolate=False, fill_bad_pix=True, quadratic_refinement=True, integrate=False, clip=None, filter=None, verbose=True, guess=False)[source]

Update the pixel_map by minimising an error metric within the search_window.

Parameters
  • data (ndarray, float32, (N, M, L)) –

    Input diffraction data \(I^{z_1}_\phi\), the \(^{z_1}\) indicates the distance between the virtual source of light and the \(_\phi\) indicates the phase of the wavefront incident on the sample surface. The dimensions are given by:

    • N = number of frames

    • M = number of pixels along the slow scan axis of the detector

    • L = number of pixels along the fast scan axis of the detector

  • mask (ndarray, bool, (M, L)) – Detector good / bad pixel mask \(M\), where True indicates a good pixel and False a bad pixel.

  • W (ndarray, float, (M, L)) – The whitefield image \(W\). This is the image one obtains without a sample in place.

  • pixel_map (ndarray, float, (2, M, L)) –

    An array containing the pixel mapping between a detector frame and the object \(ij_\text{map}\), such that:

    \[\begin{split}I^{z_1}_{\phi}[n, i, j] = W[i, j] I^\infty[&\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0,\\ &\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0]\end{split}\]

  • n0 (float) –

    Slow scan offset to the pixel mapping such that:

    \[\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0 \ge -0.5 \quad\text{for all } i,j\]

  • m0 (float) –

    Fast scan offset to the pixel mapping such that:

    \[\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0 \ge -0.5 \quad\text{for all } i,j\]

  • dij_n (ndarray, float, (N, 2)) – An array containing the sample shifts for each detector image in pixel units \(\Delta ij_n\).

  • search_window (int, len 2 sequence, optional) – The pixel mapping will be updated in a square area of side length “search_window”. If “search_window” is a length 2 sequence (e.g. [8,12]) then the search area will be rectangular with [ss_range, fs_range]. This value/s are in pixel units.

  • subpixel (bool, optional) – If True then bilinear interpolation is used to evaluate subpixel locations.

  • filter (None or float, optional) – If float then apply a gaussian filter to the pixel_maps, ignoring masked pixels. The “filter” is equal to the sigma of the Gaussian in pixel units.

  • verbose (bool, optional) – print what I’m doing.

Returns

  • pixel_map (ndarray, float, (2, M, L)) – An array containing the updated pixel mapping.

  • res (dictionary) –

    A dictionary containing diagnostic information:

    • error_mapndarray, float, (M, L)

      The minimum value of the error metric at each detector pixel

Notes

The following error metric is minimised with respect to \(\text{ij}_\text{map}[0, i, j]\):

\[\begin{split}\begin{align} \varepsilon[i, j] = \sum_n \bigg(I^{z_1}_{\phi}[n, i, j] - W[i, j] I^\infty[&\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0,\\ &\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0]\bigg)^2 \end{align}\end{split}\]
speckle_tracking.calc_error(data, mask, W, dij_n, I, pixel_map, n0, m0, subpixel=False, verbose=True)[source]
Parameters
  • data (ndarray) – Input data, of shape (N, M, L).

  • mask (ndarray) – Boolean array of shape (M, L), where True indicates a good pixel and False a bad pixel.

  • W (ndarray) – Float array of shape (M, L) containing the estimated whitefield.

  • dij_n (ndarray) – Float array of shape (N, 2) containing the object translations that have been mapped onto the detector’s frame of reference.

  • I (ndarray) – Float array of shape (U, V), this is essentially an object map.

  • pixel_map (ndarray, (2, M, L)) –

    An array containing the pixel mapping between a detector frame and the object, such that:

    \[\begin{split}I^{z_1}_{\phi}[n, i, j] = W[i, j] I^\infty[&\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0,\\ &\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0]\end{split}\]

  • n0 (float) –

    Slow scan offset to the pixel mapping such that:

    \[\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0 \ge -0.5 \quad\text{for all } i,j\]

  • m0 (float) –

    Fast scan offset to the pixel mapping such that:

    \[\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0 \ge -0.5 \quad\text{for all } i,j\]

  • subpixel (bool, optional) – If True then use bilinear subpixel interpolation non-integer pixel mappings.

  • verbose (bool, optional) – print what I’m doing.

  • minimum_overlap (float or None, optional) – Default is None. If float then the the object will be set to -1 where the number of data points contributing to that value is less than “minimum_overlap”.

Returns

  • error_total (float) – The global error value, \(\varepsilon = \sum_{n,i,j} \varepsilon[n, i, j]\).

  • error_frame (ndarray) – Float array of shape (N,). The average pixel error per detector frame, \(\varepsilon_\text{frame}[n] = \langle \varepsilon[n, i, j] \rangle_{i,j}\).

  • error_pixel (ndarray) – Float array of shape (M, L). The average pixel error per detector pixel, \(\varepsilon_\text{pixel}[i, j] = \langle \varepsilon[n, i, j]\rangle_{n}\).

Notes

The error, per pixel and per frame, is given by:

\[\begin{split}\begin{align} \varepsilon[n, i, j] = M[i,j] \bigg[ I_\Phi[n, i, j] - W[i, j] I_0[&\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0,\\ &\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0]\bigg]^2 \end{align}\end{split}\]
speckle_tracking.integrate_pixel_map(pixel_map, weight, wavelength, z, x_pixel_size, y_pixel_size, dx, dy, remove_astigmatism=False, maxiter=3000)[source]

Integrate the pixel map to produce the phase and angle profile of the lens pupil.

Notes

The pixel map (\(u[i, j]\)) defines the pixel mapping between the reference and recorded images. In 1D this is given by:

\(I[i] = O[u[i] + n_0]\)

In physical units this represents the mapping:

\[\begin{split}\begin{align} I(x) &= O(x - \nabla \Phi'(x)) \\ I[i] &= I(i \Delta_{ss}) \\ &= O(i \Delta_{ss} - \nabla \Phi'(i \Delta_{ss})) \\ &= O[i \frac{\Delta_{ss}}{\Delta_x} - \frac{1}{\Delta_x}\nabla \Phi'(i \Delta_{ss})] \; \text{where}\\ O[i] &\equiv O(i \Delta_x) \quad I[i] \equiv I(i \Delta_{ss}) \; \text{and} \\ \Phi'(x) &\equiv \frac{\lambda z}{2\pi} \Phi(x) \end{align}\end{split}\]

Therefore we have that:

\[\begin{split}\begin{align} u[i] &= i \frac{\Delta_{ss}}{\Delta_x} - \frac{1}{\Delta_x}\nabla \Phi'(i \Delta_{ss}) - n_0 \; \text{or}\\ \nabla \Phi'(i \Delta_{ss}) &= i \Delta_{ss} - \Delta_x(u[i] + n_0) \\ \Phi'(i \Delta_{ss}) &\approx \Delta_{ss} \sum_{j=0}^{j=i} \left[ j \Delta_{ss} - \Delta_x (u[j] + n_0)\right] \\ &= \Delta_x\Delta_{ss}(i+1) \left[\frac{i\Delta_{ss}}{2\Delta_x} - n_0 - \frac{1}{i+1} \sum_{j=0}^{j=i} u[j] \right] \end{align}\end{split}\]

The angle of each ray pointing from pixel i in the detector is given by \(\theta(i\Delta_{ss}) = -\frac{1}{z}\nabla \Phi'(i\Delta_{ss}) = \frac{1}{z} \left(\Delta_x(u[i] + n_0) - i \Delta_{ss}\right)\). Another way to think of this is:

\[\begin{split}\begin{align} \theta(i\Delta_{ss}) &= \frac{\text{position in object}-\text{position in detector}}{z} \\ &= \frac{(u[i] + n_0)\Delta_x - i \Delta_{ss}}{z} \\ \end{align}\end{split}\]

The residual angles, the ray angles after subtracting the global curvature \(\theta_r(x) = \theta(x) - \theta_0(x)\), are given by:

\[\begin{split}\begin{align} \nabla \Phi_0'(x) &= \frac{\lambda z}{2\pi} \frac{2\pi x}{\lambda z_r} = x \frac{z}{z_r} \\ \nabla \Phi_r'(i \Delta_{ss}) &= i \Delta_{ss} - \Delta_x(u[i] + n_0) - i \Delta_{ss} \frac{z}{z_r} \\ &= i \Delta_{ss} \frac{z_r-z}{z_r} - \Delta_x(u[i] + n_0) \\ \theta_r(x) &= \frac{\Delta_x(u[i] + n_0) - i \Delta_{ss} \frac{z_r - z}{z_r}}{z} \\ \Phi_r'(i \Delta_{ss}) &= \Delta_{ss} \sum_{j=0}^{j=i} \left[ j\frac{z_r - z}{z_r} \Delta_{ss} - \Delta_x (u[j] + n_0)\right] \end{align}\end{split}\]
speckle_tracking.docstring_glossary()[source]
Parameters
  • data (ndarray, float32, (N, M, L)) –

    Input diffraction data \(I^{z_1}_\phi\), the \(^{z_1}\) indicates the distance between the virtual source of light and the \(_\phi\) indicates the phase of the wavefront incident on the sample surface. The dimensions are given by:

    • N = number of frames

    • M = number of pixels along the slow scan axis of the detector

    • L = number of pixels along the fast scan axis of the detector

  • mask (ndarray, bool, (M, L)) – Detector good / bad pixel mask \(M\), where True indicates a good pixel and False a bad pixel.

  • x_pixel_size (float) – The side length of a detector pixel in metres \(\Delta_{ss}\), along the slow scan axis.

  • y_pixel_size (float) – The side length of a detector pixel in metres \(\Delta_{fs}\), along the fast scan axis.

  • z (float) – The distance between the focus and the detector in metres \(z = z_1 + z_2\).

  • wav (float) – The wavelength of the imaging radiation in metres \(\lambda\).

  • W (ndarray, float, (M, L)) – The whitefield image \(W\). This is the image one obtains without a sample in place.

  • roi (array_like, int, (4,)) – Length 4 list of integers e.g. roi = [10, 400, 23, 500], indicates that most of the interesting data in a frame will be in the region: frame[roi[0]:roi[1], roi[2]:roi[3]]

  • defocus (float) – The average focus to sample distance \(z_1\).

  • dz (float) –

    The difference between the average defocus and the defocus along the fast and slow axes of the detector:

    \[\begin{split}z_{ss} &= z_1 + dz \\ z_{fs} &= z_1 - dz\end{split}\]

    \(z_1 > 0\) indicates that the sample is downstream of focus.

  • O (ndarray, float, (U, V)) – This is essentially a defocused image of the object \(O\) or \(I^\infty\). It is the image one would obtain with plane wave illumination with the detector placed some distance from the sample.

  • n0 (float) –

    Slow scan offset to the pixel mapping such that:

    \[\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0 \ge -0.5 \quad\text{for all } i,j\]

  • m0 (float) –

    Fast scan offset to the pixel mapping such that:

    \[\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0 \ge -0.5 \quad\text{for all } i,j\]

  • pixel_map (ndarray, float, (2, M, L)) –

    An array containing the pixel mapping between a detector frame and the object \(ij_\text{map}\), such that:

    \[\begin{split}I^{z_1}_{\phi}[n, i, j] = W[i, j] I^\infty[&\text{ij}_\text{map}[0, i, j] - \Delta ij[n, 0] + n_0,\\ &\text{ij}_\text{map}[1, i, j] - \Delta ij[n, 1] + m_0]\end{split}\]

  • dij_n (ndarray, float, (N, 2)) – An array containing the sample shifts for each detector image in pixel units \(\Delta ij_n\).