pypret.retrieval package

This module provides classes to calculate parametrized nonlinear process spectra (PNPS), such as frequency-resolved optical gating (FROG), interferometric FROG (iFROG), dispersion scan (d-scan), time-domain ptychography (TDP) and pulse-shaper assisted methods such as multiphoton intrapulse interference phase scan (MIIPS).

The code follows the notation used in [Geib2019] and its supplement.

Currently only the abovementioned methods are implemented. But the code is written in such way that including new pulse measurement methods is very easy. If it is a method using a collinear nonlinearity, subclass from CollinearPNPS, otherwise from NoncollinearPNPS.

In the collinear case only self.mask(parameter) has to be implemented which calculates the used linear parametrization operator. In the non-collinear case the function _calculate has to be implemented which calculates and returns the PNPS trace T_mn and the PNPS signal S_mk.

Retrieval algorithms

pypret.retrieval.retriever.Retriever(pnps: pypret.pnps.BasePNPS, method: str = 'copra', maxiter=300, maxfev=None, logging=False, verbose=False, **kwargs) → pypret.retrieval.retriever.BaseRetriever[source]

Creates a retriever instance.

Parameters:
  • pnps (PNPS) – A PNPS instance that is used to simulate a PNPS measurement.
  • method (str, optional) –
    Type of solver. Should be one of
    • ’copra’ (see here)
    • ’gpa’ (see here)
    • ’gp-dscan’ (see here)
    • ’pcgpa’ (see here)
    • ’pie’ (see here)
    • ’lm’ (see here)
    • ’bfgs’ (see here)
    • ’de’ (see here)
    • ’nelder-mead’ (see here)

    ’copra’ is the default choice.

  • maxiter (int, optional) – The maximum number of algorithm iterations. The default is 300.
  • maxfev (int, optional) – The maximum number of function evaluations. If given, the algorithms stop before this number is reached. Not all algorithms support this feature. Default is None, in which case it is ignored.
  • logging (bool, optional) – Stores trace errors and pulses over the iterations if supported by the retriever class. Default is False.
  • verbose (bool, optional) – Prints out trace errors during the iteration if supported by the retriever class. Default is False.
class pypret.retrieval.step_retriever.COPRARetriever(pnps, alpha=0.25, **kwargs)[source]

This module implements the common pulse retrieval algorithm [Geib2019].

__init__(pnps, alpha=0.25, **kwargs)[source]

For a full documentation of the arguments see Retriever.

Parameters:alpha (float, optional) – Scales the step size in the global stage of COPRA. Higher values mean potentially faster convergence but less accuracy. Lower values provide higher accuracy for the cost of speed. Default is 0.25.
method = 'copra'
class pypret.retrieval.step_retriever.PCGPARetriever(pnps, decomposition='power', **kwargs)[source]

This class implements the principal components generalized projections algorithm (PCGPA) for SHG-FROG.

We follow the algorithm as described in [Kane1999] but use the PNPS formalism from [Geib2019] and some minor modifications:

  • it supports both the singular value decomposition and the power method to find/approximate the largest eigenvector.
  • the projection includes the scaling factor µ. This makes the method robust against initial guesses with the wrong magnitude. It should have no adverse effect.
__init__(pnps, decomposition='power', **kwargs)[source]

For a full documentation of the arguments see Retriever.

Parameters:decomposition (str, optional) – It specifies how the FROG signal is decomposed. If power (the default) the power method is used to find the largest eigenvalue. If svd a full singular value decomposition is performed. This is potentially much slower but more accurate.
method = 'pcgpa'
supported_schemes = ['shg-frog']
class pypret.retrieval.step_retriever.GPARetriever(pnps, step_size='exact', **kwargs)[source]

Implements the classical generalized projections algorithm for SHG-FROG as described in [DeLong1994] and [Trebino2000].

As far as I know the determination of the step size in GPA is not made explicit in the publications. It is usually done in a line search. In this implementation we offer three different options:

  • an exact line search using a Brent style minimizer
  • a backtracking (inexact) line search using the Armijo-Goldstein condition with c=0.5 and tau=0.5.
  • the same heuristic choice for the step size used in copra.

The last method is the fastest, but as the first is the classic choice for GPA, it is the default.

__init__(pnps, step_size='exact', **kwargs)[source]

For a full documentation of the arguments see Retriever.

Parameters:step_size (str, optional) – Specifies how the step size of the gradient step in GPA is determined. Default is exact which performs an exact line search. inexact performs a backtracking line search and copra uses the ad-hoc estimates for the step size used in COPRA.
method = 'gpa'
supported_schemes = ['shg-frog']
class pypret.retrieval.step_retriever.PIERetriever(pnps, logging=False, verbose=False, **kwargs)[source]

This class implements a pulse retrieval algorithm for SHG-FROG based on the ptychographical iterative engine (PIE). It is based on the paper [Sidorenko2016] and its erratum [Sidorenko2017].

We modified the algorithm to include the scaling factor µ in the projection. This makes the method robust against initial guesses with the wrong magnitude. It should have no adverse effect.

method = 'pie'
supported_schemes = ['shg-frog']
class pypret.retrieval.step_retriever.GPDSCANRetriever(pnps, logging=False, verbose=False, **kwargs)[source]

This class implements a pulse retrieval algorithm for SHG and THG d-scan based on the paper [Miranda2017].

In our tests we found that it does not converge in the noiseless case. In other words the global solution to the least-squares problem is not a fixed point of the iteration.

method = 'gp-dscan'
supported_schemes = ['shg-dscan', 'thg-dscan']
class pypret.retrieval.nlo_retriever.LMRetriever(pnps, ftol=1e-08, xtol=1e-08, gtol=1e-08, lm_verbose=0, **kwargs)[source]

Implements pulse retrieval based on the Levenberg-Marquadt algorithm.

This is an efficient nonlinear least-squares solver, however, it will still be very slow for large pulses (N > 256). The reason is that the (MN x N) Jacobian is evaluated using numerical differentiation.

The recommendation is to use this method either on small problems or to refine or verify solutions provided by a different algorithm.

__init__(pnps, ftol=1e-08, xtol=1e-08, gtol=1e-08, lm_verbose=0, **kwargs)[source]

For a full documentation of the arguments see Retriever.

For the documentation of ftol, xtol, gtol see the documentation of scipy.optimize.least_squares(). They are passed directly to the optimizer. If you want to run the optimizer for a fixed number of iterations, set all values to 1e-14 to effectively disable the stopping criteria.

method = 'lm'
class pypret.retrieval.nlo_retriever.NMRetriever(pnps, logging=False, verbose=False, **kwargs)[source]

This retriever uses the gradient-free Nelder-Mead algorithm.

method = 'nm'
class pypret.retrieval.nlo_retriever.DERetriever(pnps, logging=False, verbose=False, **kwargs)[source]

This retriever uses the gradient-free differential evolution algorithm.

It tries to match the parameters described in [Escoto2018] as far as they are mentioned. No further effort was made to optimize them. If you are interested in using DE as a pulse retrieval algorithm you are advised to study the documentation at scipy.optimize.differential_evolution().

The initial population in our implementation is based on the provided guess with added complex, Gaussian noise of 5% of the maximum amplitude. In our tests we saw no convergence when starting from completely random initial guesses.

method = 'de'
class pypret.retrieval.nlo_retriever.BFGSRetriever(pnps, logging=False, verbose=False, **kwargs)[source]

This retriever uses the BFGS algorithm with numerical differentiation.

method = 'bfgs'

API

class pypret.retrieval.retriever.BaseRetriever(pnps, logging=False, verbose=False, **kwargs)[source]

The abstract base class for pulse retrieval.

This class implements common functionality for different retrieval algorithms.

__init__(pnps, logging=False, verbose=False, **kwargs)[source]

Initialize self. See help(type(self)) for accurate signature.

method = None
result(pulse_original=None, full=True)[source]

Analyzes the retrieval results in one retrieval instance and processes it for plotting or storage.

retrieve(measurement, initial_guess, weights=None, **kwargs)[source]

Retrieve pulse from measurement starting at initial_guess.

Parameters:
  • measurement (MeshData) – A MeshData instance that contains the PNPS measurement. The first axis has to correspond to the PNPS parameter, the second to the frequency. The data has to be the measured _intensity_ over the frequency (not wavelength!). The second axis has to match exactly the frequency axis of the underlying PNPS instance. No interpolation is done.
  • initial_guess (1d-array) – The spectrum of the pulse that is used as initial guess in the iterative retrieval.
  • weights (1d-array) – Weights that are attributed to the measurement for retrieval. In the case of (assumed) Gaussian uncertainties with standard deviation sigma they should correspond to 1/sigma. Not all algorithms support using the weights.
  • kwargs (dict) – Can override retrieval options specified in __init__().

Notes

This function provides no interpolation or data processing. You have to write a retriever wrapper for that purpose.

supported_schemes = None
trace_error(spectrum, store=True)[source]

Calculates the trace error from the pulse spectrum.

class pypret.retrieval.step_retriever.StepRetriever(pnps, logging=False, verbose=False, **kwargs)[source]