solvers#

Contains different solver implementations.

class tike.ptycho.solvers.DmOptions[source]#

Bases: IterativeOptions

DmOptions(num_batch: ‘int’ = 1, batch_method: ‘str’ = ‘wobbly_center’, rescale_method: ‘str’ = ‘mean_of_abs_object’, rescale_period: ‘int’ = 10, num_iter: ‘int’ = 1, convergence_window: ‘int’ = 0, time_limit: ‘float’ = numpy.inf)

__init__(num_batch=1, batch_method='wobbly_center', rescale_method='mean_of_abs_object', rescale_period=10, num_iter=1, convergence_window=0, time_limit=numpy.inf)#
Parameters
  • num_batch (int) –

  • batch_method (str) –

  • rescale_method (str) –

  • rescale_period (int) –

  • num_iter (int) –

  • convergence_window (int) –

  • time_limit (float) –

Return type

None

batch_method: str = 'wobbly_center'#

The name of the batch selection method. Choose from the cluster methods in the tike.cluster module.

convergence_window: int = 0#

The number of epochs to consider for convergence monitoring. Set to any value less than 2 to disable.

costs: List[List[float]]#

The objective function value at previous iterations. One list is returned for each mini-batch.

name: str = 'dm'#

The name of the algorithm.

num_batch: int = 1#

The dataset is divided into this number of groups where each group is processed simultaneously.

num_iter: int = 1#

The number of epochs to process before returning.

rescale_method: str = 'mean_of_abs_object'#

How we control object/probe scaling in the ptycho optimization problem.

The options here are:

‘mean_of_abs_object’ = The default is using the constraint that the mean of the absolute value of the object must be approx 1.0, which will then rescale the object and probe so that this constraint is true.

‘constant_probe_photons’ = Rescale the shared probe modes so that its intensity L2 norm equals to some number of photons specified elsewhere (e.g. toml file). If not specified elsewhere, the average (wrt scan positions) L2 norm of the diffraction intensity measurements is used.

rescale_period: int = 10#

How often we control object/probe scaling in the ptycho optimization problem.

The default is perform rescaling of the object/probe every 10 epochs

times: List[float]#

The per-iteration wall-time for each previous iteration.

class tike.ptycho.solvers.LstsqOptions[source]#

Bases: IterativeOptions

LstsqOptions(num_batch: ‘int’ = 1, batch_method: ‘str’ = ‘wobbly_center’, rescale_method: ‘str’ = ‘mean_of_abs_object’, rescale_period: ‘int’ = 10, num_iter: ‘int’ = 1, convergence_window: ‘int’ = 0, time_limit: ‘float’ = numpy.inf)

__init__(num_batch=1, batch_method='wobbly_center', rescale_method='mean_of_abs_object', rescale_period=10, num_iter=1, convergence_window=0, time_limit=numpy.inf)#
Parameters
  • num_batch (int) –

  • batch_method (str) –

  • rescale_method (str) –

  • rescale_period (int) –

  • num_iter (int) –

  • convergence_window (int) –

  • time_limit (float) –

Return type

None

batch_method: str = 'wobbly_center'#

The name of the batch selection method. Choose from the cluster methods in the tike.cluster module.

convergence_window: int = 0#

The number of epochs to consider for convergence monitoring. Set to any value less than 2 to disable.

costs: List[List[float]]#

The objective function value at previous iterations. One list is returned for each mini-batch.

name: str = 'lstsq_grad'#

The name of the algorithm.

num_batch: int = 1#

The dataset is divided into this number of groups where each group is processed sequentially.

num_iter: int = 1#

The number of epochs to process before returning.

rescale_method: str = 'mean_of_abs_object'#

How we control object/probe scaling in the ptycho optimization problem.

The options here are:

‘mean_of_abs_object’ = The default is using the constraint that the mean of the absolute value of the object must be approx 1.0, which will then rescale the object and probe so that this constraint is true.

‘constant_probe_photons’ = Rescale the shared probe modes so that its intensity L2 norm equals to some number of photons specified elsewhere (e.g. toml file). If not specified elsewhere, the average (wrt scan positions) L2 norm of the diffraction intensity measurements is used.

rescale_period: int = 10#

How often we control object/probe scaling in the ptycho optimization problem.

The default is perform rescaling of the object/probe every 10 epochs

times: List[float]#

The per-iteration wall-time for each previous iteration.

class tike.ptycho.solvers.PtychoParameters[source]#

Bases: object

A class for storing the ptychography forward model parameters.

New in version 0.22.0.

__init__(probe, psi, scan, eigen_probe=None, eigen_weights=None, algorithm_options=<factory>, exitwave_options=None, probe_options=None, object_options=None, position_options=None)#
Parameters
  • probe (numpy.typing.NDArray.numpy.csingle) –

  • psi (numpy.typing.NDArray.numpy.csingle) –

  • scan (numpy.typing.NDArray.numpy.single) –

  • eigen_probe (Optional[numpy.typing.NDArray.numpy.csingle]) –

  • eigen_weights (Optional[numpy.typing.NDArray.numpy.single]) –

  • algorithm_options (IterativeOptions) –

  • exitwave_options (Optional[ExitWaveOptions]) –

  • probe_options (Optional[ProbeOptions]) –

  • object_options (Optional[ObjectOptions]) –

  • position_options (Optional[PositionOptions]) –

Return type

None

algorithm_options: IterativeOptions#

A class containing algorithm specific parameters

eigen_probe: Optional[numpy.typing.NDArray.numpy.csingle] = None#

(EIGEN, SHARED, WIDE, HIGH) complex64 The eigen probes for all positions.

eigen_weights: Optional[numpy.typing.NDArray.numpy.single] = None#

(POSI, EIGEN, SHARED) float32 The relative intensity of the eigen probes at each position.

exitwave_options: Optional[ExitWaveOptions] = None#

A class containing settings related to exitwave updates.

object_options: Optional[ObjectOptions] = None#

A class containing settings related to object updates.

position_options: Optional[PositionOptions] = None#

A class containing settings related to position correction.

probe: numpy.typing.NDArray.numpy.csingle#

(1, 1, SHARED, WIDE, HIGH) complex64 The shared illumination function amongst all positions.

probe_options: Optional[ProbeOptions] = None#

A class containing settings related to probe updates.

psi: numpy.typing.NDArray.numpy.csingle#

(WIDE, HIGH) complex64 The wavefront modulation coefficients of the object.

resample(factor, interp)[source]#

Return a new PtychoParameter with the parameters rescaled.

Parameters
  • factor (float) –

  • interp (Optional[Callable[[numpy.ndarray, float], numpy.array]]) –

Return type

PtychoParameters

scan: numpy.typing.NDArray.numpy.single#

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

class tike.ptycho.solvers.RpieOptions[source]#

Bases: IterativeOptions

RpieOptions(num_batch: ‘int’ = 5, batch_method: ‘str’ = ‘wobbly_center’, rescale_method: ‘str’ = ‘mean_of_abs_object’, rescale_period: ‘int’ = 10, num_iter: ‘int’ = 1, convergence_window: ‘int’ = 0, time_limit: ‘float’ = numpy.inf, alpha: ‘float’ = 0.05)

__init__(num_batch=5, batch_method='wobbly_center', rescale_method='mean_of_abs_object', rescale_period=10, num_iter=1, convergence_window=0, time_limit=numpy.inf, alpha=0.05)#
Parameters
  • num_batch (int) –

  • batch_method (str) –

  • rescale_method (str) –

  • rescale_period (int) –

  • num_iter (int) –

  • convergence_window (int) –

  • time_limit (float) –

  • alpha (float) –

Return type

None

alpha: float = 0.05#

A hyper-parameter which controls the step length. RPIE becomes EPIE when this parameter is 1.

batch_method: str = 'wobbly_center'#

The name of the batch selection method. Choose from the cluster methods in the tike.cluster module.

convergence_window: int = 0#

The number of epochs to consider for convergence monitoring. Set to any value less than 2 to disable.

costs: List[List[float]]#

The objective function value at previous iterations. One list is returned for each mini-batch.

name: str = 'rpie'#

The name of the algorithm.

num_batch: int = 5#

The dataset is divided into this number of groups where each group is processed sequentially.

num_iter: int = 1#

The number of epochs to process before returning.

rescale_method: str = 'mean_of_abs_object'#

How we control object/probe scaling in the ptycho optimization problem.

The options here are:

‘mean_of_abs_object’ = The default is using the constraint that the mean of the absolute value of the object must be approx 1.0, which will then rescale the object and probe so that this constraint is true.

‘constant_probe_photons’ = Rescale the shared probe modes so that its intensity L2 norm equals to some number of photons specified elsewhere (e.g. toml file). If not specified elsewhere, the average (wrt scan positions) L2 norm of the diffraction intensity measurements is used.

rescale_period: int = 10#

How often we control object/probe scaling in the ptycho optimization problem.

The default is perform rescaling of the object/probe every 10 epochs

times: List[float]#

The per-iteration wall-time for each previous iteration.

tike.ptycho.solvers.crop_fourier_space(x, w)[source]#

Crop x assuming a 2D frequency space image with zero frequency in corner.

Parameters
  • x (numpy.ndarray) –

  • w (int) –

Return type

numpy.ndarray

tike.ptycho.solvers.dm(op, comm, data, batches, *, parameters, epoch)[source]#

Solve the ptychography problem using the difference map approach.

Parameters
  • op (tike.operators.Ptycho) – A ptychography operator.

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

  • data (list((FRAME, WIDE, HIGH) float32, ...)) – A list of unique CuPy arrays for each device containing the intensity (square of the absolute value) of the propagated wavefront; i.e. what the detector records. FFT-shifted so the diffraction peak is at the corners.

  • batches (list(list((BATCH_SIZE, ) int, ...), ...)) – A list of list of indices along the FRAME axis of data for each device which define the batches of data to process simultaneously.

  • parameters (tike.ptycho.solvers.PtychoParameters) – An object which contains reconstruction parameters.

  • epoch (int) –

Returns

result – An object which contains reconstruction parameters.

Return type

tike.ptycho.solvers.PtychoParameters

References

Thibault, Pierre, Martin Dierolf, Oliver Bunk, Andreas Menzel, and Franz Pfeiffer. “Probe retrieval in ptychographic coherent diffractive imaging.” Ultramicroscopy 109, no. 4 (2009): 338-343.

See also

tike.ptycho

tike.ptycho.solvers.lstsq_grad(op, comm, data, batches, *, parameters, epoch)[source]#

Solve the ptychography problem using Odstrcil et al’s approach.

Object and probe are updated simultaneously using optimal step sizes computed using a least squares approach.

Parameters
  • op (tike.operators.Ptycho) – A ptychography operator.

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

  • data (list((FRAME, WIDE, HIGH) float32, ...)) – A list of unique CuPy arrays for each device containing the intensity (square of the absolute value) of the propagated wavefront; i.e. what the detector records. FFT-shifted so the diffraction peak is at the corners.

  • batches (list(list((BATCH_SIZE, ) int, ...), ...)) – A list of list of indices along the FRAME axis of data for each device which define the batches of data to process simultaneously.

  • parameters (tike.ptycho.solvers.PtychoParameters) – An object which contains reconstruction parameters.

  • epoch (int) –

Returns

result – A dictionary containing the updated keyword-only arguments passed to this function.

Return type

dict

References

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

See also

tike.ptycho

tike.ptycho.solvers.rpie(op, comm, data, batches, *, parameters, epoch)[source]#

Solve the ptychography problem using regularized ptychographical engine.

The rPIE update direction can be shown to be equivalent to a conventional gradient descent direction but rescaled by the preconditioner term. i.e. If the rPIE step size (alpha) is 0 and the preconditioner is zero, we have the vanilla gradient descent direction.

Parameters
  • op (tike.operators.Ptycho) – A ptychography operator.

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

  • data (list((FRAME, WIDE, HIGH) float32, ...)) – A list of unique CuPy arrays for each device containing the intensity (square of the absolute value) of the propagated wavefront; i.e. what the detector records. FFT-shifted so the diffraction peak is at the corners.

  • batches (list(list((BATCH_SIZE, ) int, ...), ...)) – A list of list of indices along the FRAME axis of data for each device which define the batches of data to process simultaneously.

  • parameters (tike.ptycho.solvers.PtychoParameters) – An object which contains reconstruction parameters.

  • epoch (int) –

Returns

result – An object which contains the updated reconstruction parameters.

Return type

tike.ptycho.solvers.PtychoParameters

References

Maiden, Andrew M., and John M. Rodenburg. 2009. “An Improved Ptychographical Phase Retrieval Algorithm for Diffractive Imaging.” Ultramicroscopy 109 (10): 1256–62. https://doi.org/10.1016/j.ultramic.2009.05.012.

Andrew Maiden, Daniel Johnson, and Peng Li, “Further improvements to the ptychographical iterative engine,” Optica 4, 736-745 (2017) https://doi.org/10.1364/OPTICA.4.000736

See also

tike.ptycho

tike.ptycho.solvers.update_preconditioners(comm, operator, scan, probe, psi, object_options=None, probe_options=None)[source]#

Update the probe and object preconditioners.

Parameters
Return type

Tuple[ObjectOptions, ProbeOptions]