Gradient-based Engines

Gradient-based Engines#

Learning more about gradient-based reconstruction engines/algorithms such as the maximum likelihood algorithm.

In contrast to the projectional and stochastic algorithms, the maximum-likelihood approach allows to explicitely define a noise model for the ptychographic reconstruction. For ptychography, the negative log-likelihood can be written as

\[ \mathcal{L} = - \log \prod_k \prod_q p(I_{kq}|PO) \]

where \(p(I_{kq}|PO)\) is the probability of measuring a photon \(I_{kq}\) at scan view \(k\) and detector pixel \(q\) given probe \(P\) and object \(O\). For a Gaussian distribution the likelihood reduces to a sum of least squares

\[ \mathcal{L} = \sum_k \sum_q \frac{(|\Psi_{kq}|^2 - I_{kq})^2}{2\sigma^2_{kq}} \]

which can be used to minimise with respect to \(P\) and \(O\) using a non-linear conjugate gradient (CG) approach by defining the gradient

\[ g = \left(\frac{\partial \mathcal{L}}{\partial P}, \frac{\partial\mathcal{L}}{\partial O}\right) \]

which in turn can be used to define the new conjugate search direction for the \(n-th\) iteration step

\[ \Delta^{(n)} = -g^{(n)} + \beta^{(n)}\Delta^{(n-1)} \]

where \(\beta^{(n)}\) can be calculated via the Polak-Ribiere formula. For more details about this maximum-likelihood based CG algorithm, see Thibault P. and Guizar-Sicairos M. (2012) which also includes a description of preconditioners and regularisers that can help with the convergence of the algorithm. In particular, the scale preconditioner can be useful for weakly scattering objects - assuming that the object deviates only weakly from unity - by enforcing a near equality between the probe and object gradients. This parameter can be controlled using scale_precond which is True by default and scale_probe_object which is \(1.0\) by default.

p.engines.engine00.scale_precond = True
p.engines.engine00.scale_probe_object = 1.

In addition, regularisers can be used to reduce difficulties that can arise from ill-posed problems such as the inverse phase problem in ptychography. To impose some form of smoothness in the final reconstruction, we can use a simple Gaussian regulariser which can be toggled using reg_del2 (True by default) and tuned using reg_del2_amplitude which defaults to \(0.01\). In this moonflower example, we are changing the amplitude to \(1.\) for a more aggressive smoothing strategy.

p.engines.engine00.reg_del2 = True 
p.engines.engine00.reg_del2_amplitude = 1.

Challenge

Modify different ML engine parameters and observe their impact on the outcome of the reconstruction.


import ptypy
import ptypy.utils as u

p = u.Param()
p.verbose_level = "interactive"
p.io = u.Param()
p.io.rfile = None
p.io.autosave = u.Param(active=False)
p.io.interaction = u.Param(active=False)

# Live-plotting
p.io.autoplot = u.Param()
p.io.autoplot.active=True
p.io.autoplot.threaded = False
p.io.autoplot.layout = "jupyter"
p.io.autoplot.interval = 10

p.scans = u.Param()
p.scans.MF = u.Param()
p.scans.MF.name = "Full"
p.scans.MF.data= u.Param()
p.scans.MF.data.name = "MoonFlowerScan"
p.scans.MF.data.shape = 128
p.scans.MF.data.num_frames = 200
p.scans.MF.data.save = None
p.scans.MF.data.density = 0.2
p.scans.MF.data.photons = 1e8
p.scans.MF.data.psf = 0.

# Define reconstruction engines
p.engines = u.Param()

# Maximum Likelihood (ML) engine
p.engines.engine00 = u.Param()
p.engines.engine00.name = "ML"
p.engines.engine00.ML_type = "Gaussian"
p.engines.engine00.numiter = 300
p.engines.engine00.numiter_contiguous = 10
p.engines.engine00.reg_del2 = True 
p.engines.engine00.reg_del2_amplitude = 1.
p.engines.engine00.scale_precond = True
p.engines.engine00.scale_probe_object = 1.

P = ptypy.core.Ptycho(p,level=5)