probe#

Functions related to creating and manipulating probe arrays.

Ptychographic probes are represented as two separate components: a shared probe whose values are the same for all positions and the varying component. The former is required as it provides the shared probe constraint for ptychography and the later relaxes the former constraint to accomodate real-world illuminations which may vary with time.

The shared component consist of a single array representing at least one probe each of which may have an accompanying varying component.

The varying components are stored sparsely as two arrays, and the full representation of the varying comonents are only combined as needed. The first array is an array of eigen probes (principal components) spanning the space of the probe variation of all positions and the second is an array of weights that map the variation for each position into this space.

Each probe may have its own set of eigen probes. The unique probe at a given position is reconstructed by adding the shared probe to a weighted sum of the eigen probes.

varying_probe = weights[0] * probe + np.sum(weights[1:] * eigen_probes)

Design comments#

In theory, the probe representation could be implemented in as little as two arrays: one with all of the shared components where the probe becomes the first eigen probe and and one with the weights. Choosing to keep the eigen probes separate from the probe as a third array provides backwards compatability and allows for storing fewer eigen probes in the case when only some probes are allowed to vary.

class tike.ptycho.probe.ProbeOptions[source]#

Bases: object

Manage data and setting related to probe correction.

__init__(recover_probe=False, update_start=0, update_period=1, init_rescale_from_measurements=True, probe_photons=numpy.nan, force_orthogonality=False, force_centered_intensity=False, force_sparsity=0.0, use_adaptive_moment=False, vdecay=0.999, mdecay=0.9, probe_support=0.0, probe_support_radius=0.35, probe_support_degree=2.5, additional_probe_penalty=0.0, median_filter_abs_probe=False, median_filter_abs_probe_px=(1.0, 1.0))#
Parameters
  • recover_probe (bool) –

  • update_start (int) –

  • update_period (int) –

  • init_rescale_from_measurements (bool) –

  • probe_photons (float) –

  • force_orthogonality (bool) –

  • force_centered_intensity (bool) –

  • force_sparsity (float) –

  • use_adaptive_moment (bool) –

  • vdecay (float) –

  • mdecay (float) –

  • probe_support (float) –

  • probe_support_radius (float) –

  • probe_support_degree (float) –

  • additional_probe_penalty (float) –

  • median_filter_abs_probe (bool) –

  • median_filter_abs_probe_px (Tuple[float, float]) –

Return type

None

additional_probe_penalty: float = 0.0#

Penalty applied to the last probe for existing.

This penalty encourages the probe energy to concentrate in the lower order modes. The penalty starts at zero for the first probe and increases linearly to this value. For example, for three probes, the penalties aplied are [0.0, 0.5, 1.0].

This is a soft constraint as opposed to force_sparsity which is a hard constraint.

copy_to_device(comm)[source]#

Copy to the current GPU memory.

Return type

ProbeOptions

copy_to_host()[source]#

Copy to the host CPU memory.

Return type

ProbeOptions

force_centered_intensity: bool = False#

Forces the probe intensity to be centered.

force_orthogonality: bool = False#

Forces probes to be orthogonal each iteration.

force_sparsity: float = 0.0#

Forces this proportion of zero elements.

init_rescale_from_measurements: bool = True#

Initial rescaling of probe using measured intensity.

m: Optional[numpy.typing.NDArray]#

The first moment for adaptive moment.

mdecay: float = 0.9#

The proportion of the first moment that is previous first moments.

median_filter_abs_probe: bool = False#

Binary switch on whether to apply a median filter to absolute value of each shared probe mode.

median_filter_abs_probe_px: Tuple[float, float] = (1.0, 1.0)#

A 2-element tuple with the median filter pixel widths along each dimension.

power: List[List[float]]#

The power of the primary probe modes at each iteration.

preconditioner: Optional[numpy.typing.NDArray]#
probe_photons: float#

The shared probe mode intensity must add up to this number.

If we do not give a number for this in the parameters.toml file, then it will default to the average of the measurement intensity scaling.

probe_support: float = 0.0#

Weight of the finite probe support constraint; zero or greater.

This support constraint encourages round probes energy concentrated at the center of the probe grid. Higher support increases the effect.

probe_support_degree: float = 2.5#

Degree of the supergaussian defining the probe support; zero or greater.

Controls how hard the penalty transition is outside of the radius. Degree = 0 is a flat penalty. Degree > 0, < 1 is flatter than a gaussian. Degree 1 is a gaussian. Degree > 1 is more like a top-hat than a gaussian.

probe_support_radius: float = 0.35#

Radius of finite probe support as fraction of probe grid. [0.0, 0.5].

probe_update_sum: Optional[numpy.typing.NDArray]#

Used for momentum updates.

recover_probe: bool = False#

Boolean switch used to indicate whether to update probe or not.

resample(factor, interp)[source]#

Return a new ProbeOptions with the parameters rescaled.

Parameters

factor (float) –

Return type

ProbeOptions

update_period: int = 1#

The number of epochs between probe updates

update_start: int = 0#

Start probe updates at this epoch.

use_adaptive_moment: bool = False#

Whether or not to use adaptive moment.

v: Optional[numpy.typing.NDArray]#

The second moment for adaptive moment.

vdecay: float = 0.999#

The proportion of the second moment that is previous second moments.

tike.ptycho.probe.add_modes_cartesian_hermite(probe, nmodes)[source]#

Create more probes from a 2D Cartesian Hermite basis functions.

Starting with the given probe, new modes are computed by multiplying it with a set of 2D Cartesian Hermite functions. The probes are then orthonormalized.

Parameters
  • probe ((..., 1, WIDTH, HEIGHT)) – A single probe basis.

  • nmodes (int > 0) – The number of desired probes.

Returns

probe – New probes basis.

Return type

(…, nmodes, WIDTH, HEIGHT)

References

Michal Odstrcil, Andreas Menzel, and Manuel Guizar-Sicaros. Iterative least-squares solver for generalized maximum-likelihood ptychography. Optics Express. 2018.

tike.ptycho.probe.add_modes_random_phase(probe, nmodes)[source]#

Initialize additional probe modes by phase shifting the first mode.

Parameters
  • probe ((..., M, :, :) array) – A probe with M > 0 incoherent modes.

  • nmodes (int) – The number of desired modes.

References

M. Odstrcil, P. Baksh, S. A. Boden, R. Card, J. E. Chad, J. G. Frey, W. S. Brocklesby, “Ptychographic coherent diffractive imaging with orthogonal probe relaxation.” Opt. Express 24, 8360 (2016). doi: 10.1364/OE.24.008360

tike.ptycho.probe.adjust_probe_power(probe, power=None)[source]#

Rescale the probes according to given power.

If no power is given, then probes rescaled as 1/N.

Parameters
  • probe ((..., M, :, :) array) – A probe with M > 0 incoherent modes.

  • power ((..., M, ) array) – The relative power of the probe modes.

tike.ptycho.probe.apply_median_filter_abs_probe(probe, med_filt_px)[source]#

Apply a median filter to each of the shared probe modes.

This is meant as a “quick fix” to the probe “hot spots” numerical artifacts that arise due to the sample * probe scaling ambiguity we have in phase retrieval under the projection approximation.

Parameters
  • probe (( 1, 1, SHARED, WIDE, HIGH) complex64) – The shared probes amongst all positions.

  • med_filt_px (Tuple[float, float]) – A two element array like (tuple or list) with the median filter pixel widths

Returns

probe

Return type

the filtered probe modified in-place.

tike.ptycho.probe.constrain_center_peak(probe)[source]#

Force the peak illumination intensity to the center of the probe grid.

After smoothing the intensity of the combined illumination with a gaussian filter with standard deviation sigma, the probe is shifted such that the maximum intensity is centered.

tike.ptycho.probe.constrain_probe_sparsity(probe, f)[source]#

Constrain the probe intensity so at least f fraction elements are zero.

tike.ptycho.probe.constrain_variable_probe(comm, variable_probe, weights)[source]#

Add the following constraints to variable probe weights

  1. Remove outliars from weights

  2. Enforce orthogonality once per epoch

  3. Sort the variable probes by their total energy

  4. Normalize the variable probes so the energy is contained in the weight

tike.ptycho.probe.finite_probe_support(probe, *, radius=0.5, degree=5, p=1.0)[source]#

Returns a supergaussian penalty function for finite probe support.

A mask which provides an illumination penalty is determined by the equation:

penalty = p - p * exp( -( (x / radius)**2 + (y / radius)**2 )**degree)

where the maximum penalty is p and the minium penalty is 0. This penalty function is used in the probe gradient to supress values in the probe grid far from the center. The penalty is 0 near the center and p at the edge.

Parameters
  • radius (float (0, 0.5]) – The radius of the supergaussian.

  • degree (float >= 0) – The exponent of the terms in the supergaussian equation. Controls how hard the penalty transition is outside of the radius. Degree = 0 is a flat penalty. Degree > 0, < 1 is flatter than a gaussian. Degree 1 is a gaussian. Degree > 1 is more like a top-hat than a gaussian.

tike.ptycho.probe.gaussian(size, rin=0.8, rout=1.0)[source]#

Return a complex gaussian probe distribution.

Illumination probe represented on a 2D regular grid.

A finite-extent circular shaped probe is represented as a complex wave. The intensity of the probe is maximum at the center and damps to zero at the borders of the frame.

Parameters
  • size (int) – The side length of the distribution

  • rin (float [0, 1) < rout) – The inner radius of the distribution where the dampening of the intensity will start.

  • rout (float (0, 1] > rin) – The outer radius of the distribution where the intensity will reach zero.

tike.ptycho.probe.get_varying_probe(shared_probe, eigen_probe=None, weights=None)[source]#

Construct the varying probes.

Combines shared and eigen probes with weights to return a unique probe at each scanning position.

Parameters
  • shared_probe ((..., 1, 1, SHARED, WIDE, HIGH) complex64) – The shared probes amongst all positions.

  • eigen_probe ((..., 1, EIGEN, SHARED, WIDE, HIGH) complex64) – The eigen probes for all positions.

  • weights ((..., POSI, EIGEN + 1, SHARED) float32) – The relative intensity of the eigen probes at each position.

Returns

unique_probes

Return type

(…, POSI, 1, 1, WIDE, HIGH)

tike.ptycho.probe.init_varying_probe(scan, shared_probe, num_eigen_probes, probes_with_modes=1)[source]#

Initialize arrays varying probe / eigen probes.

If num_eigen_probes is 1, then the shared probe is allowed to vary but no additional eigen probes are created.

Parameters
  • shared_probe ((..., 1, 1, SHARED, WIDE, HIGH) complex64) – The shared probes amongst all positions.

  • scan ((..., POSI, 2) float32) – The eigen probes for all positions.

  • num_eigen_probes (int) – The number of principal components used to represent the varying probe illumination.

  • probes_with_modes (int) – The number of probes that are allowed to vary.

Returns

  • eigen_probe ((…, 1, EIGEN - 1, probes_with_modes, WIDE, HIGH) complex64) – The eigen probes for all positions. None if EIGEN <= 1.

  • weights ((…, POSI, EIGEN, SHARED) float32) – The relative intensity of the eigen probes at each position. None if EIGEN < 1.

tike.ptycho.probe.orthogonalize_eig(x)[source]#

Orthogonalize modes of x using eigenvectors of the pairwise dot product.

Parameters

x ((..., nmodes, :, :) array_like complex64) – An array of the probe modes vectorized

Returns

  • x (array_like) – The orthogonalized probes

  • power (array_like) – The power of each probe

Return type

Tuple[numpy.typing.NDArray, numpy.typing.NDArray]

References

M. Odstrcil, P. Baksh, S. A. Boden, R. Card, J. E. Chad, J. G. Frey, W. S. Brocklesby, “Ptychographic coherent diffractive imaging with orthogonal probe relaxation.” Opt. Express 24, 8360 (2016). doi: 10.1364/OE.24.008360

tike.ptycho.probe.power(probe)[source]#

Return the power of each probe mode.

Parameters

probe (numpy.typing.NDArray) –

Return type

numpy.typing.NDArray

tike.ptycho.probe.rescale_probe_using_fixed_intensity_photons(probe, Nphotons, probe_power_fraction=None)[source]#

Rescale the shared probes so the sum of their intensities is Nphotons.

Parameters
  • Nphotons (float (0, inf)) – The total number of photons in the shared probe mode intensity, i.e. the sum of the intensity of the shared probe modes.

  • probe_power_fraction (array_like) – A vector of length N_p (N_p = number of shared probe modes) that contains the relative energy of each mode; must add up to 1.0

tike.ptycho.probe.simulate_varying_weights(scan, eigen_probe)[source]#

Generate weights for eigen probe that follow random sinusoid.

The amplitude of the of weights is 1, the phase shift is (0, 2π], and the period is at most one full scan.

tike.ptycho.probe.update_eigen_probe(comm, R, eigen_probe, weights, patches, diff, batches, *, batch_index, β=0.1, c=1, m=0)[source]#

Update eigen probes using residual probe updates.

This update is copied from the source code of ptychoshelves. It is similar to, but not the same as, equation (31) described by Odstrcil et al (2018). It is also different from updates described in Odstrcil et al (2016). However, they all aim to correct for probe variation.

Parameters
  • comm (tike.communicators.Comm) – An object which manages communications between both GPUs and nodes.

  • R ((POSI, 1, 1, WIDE, HIGH) complex64) – Residual probe updates; what’s left after subtracting the shared probe update from the varying probe updates for each position

  • patches ((POSI, 1, 1, WIDE, HIGH) complex64) –

  • diff ((POSI, 1, SHARED, WIDE, HIGH) complex64) –

  • eigen_probe ((1, EIGEN, SHARED, WIDE, HIGH) complex64) – The eigen probe being updated.

  • β (float) – A relaxation constant that controls how quickly the eigen probe modes are updated. Recommended to be < 1 for mini-batch updates.

  • weights ((POSI, EIGEN, SHARED) float32) – A vector whose elements are sums of the previous optimal updates for each posiiton.

References

M. Odstrcil, P. Baksh, S. A. Boden, R. Card, J. E. Chad, J. G. Frey, W. S. Brocklesby, “Ptychographic coherent diffractive imaging with orthogonal probe relaxation.” Opt. Express 24, 8360 (2016). doi: 10.1364/OE.24.008360

Michal Odstrcil, Andreas Menzel, and Manuel Guizar-Sicaros. Iterative least-squares solver for generalized maximum-likelihood ptychography. Optics Express. 2018.