SpeckleTracking#
- class pyrost.SpeckleTracking(data, dj_pix, di_pix, ds_x, ds_y, num_threads, parent, pixel_map, test_ratio, whitefield, reference_image=None, ref_orig=None, scale_map=None, test_mask=None, train_mask=None, hval=None, error=None, initial=None)#
Wrapper class for the Robust Speckle Tracking algorithm. Provides an interface to perform the reference image and lens wavefront reconstruction.
- Parameters
parent (
ReferenceType
) – The Speckle tracking data container, from which the object is derived.kwargs – Dictionary of the attributes’ data specified in attr_set and init_set.
Notes
Necessary attributes:
data : Measured frames.
di_pix : The sample’s translations along the vertical axis in pixels.
dj_pix : The sample’s translations along the horizontal axis in pixels.
ds_y : reference_image sampling interval in pixels along the vertical axis.
ds_x : reference_image sampling interval in pixels along the horizontal axis.
num_threads : Specify number of threads that are used in all the calculations.
parent : The parent
STData
container.pixel_map : The pixel mapping between the data at the detector’s plane and the reference image at the reference plane.
scale_map : Huber scaling map.
whitefield : Measured frames’ white-field.
Optional attributes:
error : Average error of the reference image and pixel mapping fit.
hval : Smoothing kernel bandwidth used in reference_image regression. The value is given in pixels.
m0 : The lower bounds of the horizontal detector axis of the reference image at the reference frame in pixels.
n0 : The lower bounds of the vertical detector axis of the reference image at the reference frame in pixels.
reference_image : The unabberated reference image of the sample.
- Raises
ValueError – If an attribute specified in attr_set has not been provided.
See also
pyrost.bin.KR_reference()
: Full details of the reference_image update using kernel regression.pyrost.bin.LOWESS_reference()
: Full details of the reference_image update using local weighted linear regression.pyrost.bin.pm_gsearch()
: Full details of the pixel_map grid search update.
- CV(hval, method='KerReg')#
Return cross-validation error for the given kernel bandwidth hval. The cross-validation error metric is calculated as follows: (i) generate a reference profile based on the training subset and (ii) find the mean-squared-error (MSE) for the test subset.
- Parameters
- Return type
- Returns
Cross-validation error.
- Raises
ValueError – If method keyword value is not valid.
See also
pyrost.bin.pm_total_error()
: Full details of the error metric.
- CV_curve(harr, method='KerReg', verbose=True)#
Return a set of cross-validation errors for a set of kernel bandwidths.
- Parameters
- Return type
- Returns
An array of cross-validation errors.
See also
pyrost.bin.pm_total_error()
: Full details of the error metric.
- create_initial()#
Create a
SpeckleTracking
object with preliminary approximation of the pixel mapping and reference profile. The object is saved into initial attribute. Necessary to calculate normalized error metrics.- Return type
- Returns
A new
SpeckleTracking
object with the preliminarySpeckleTracking
object saved in the initial attribute.
- error_profile(kind='pixel_map')#
Return a error metrics for the reference profile and pixel mapping updates. The error metrics may be normalized by the error of the initial estimations of the reference profile and pixel mapping function. The normalization is performed if
SpeckleTracking.create_initial()
was invoked before.- Parameters
kind (
str
) – Choose between generating the error metric of the pixel mapping (‘pixel_map’) or reference image (‘reference’)update.- Raises
ValueError – If kind keyword value is not valid.
- Return type
- Returns
Residual profile.
See also
SpeckleTracking.update_reference()
: More details on the reference profile update procedure.SpeckleTracking.update_pixel_map()
: More details on the pixel mapping update procedure.
- find_hopt(h0=1.0, method='KerReg', beta=0.9, epsilon=0.001, maxiter=10, gtol=1e-05, verbose=False)#
Find the optimal kernel bandwidth by finding the bandwidth, that minimized the Cross-validation error metric. The minimization process is performed with the quasi-Newton method of Broyden, Fletcher, Goldfarb, and Shanno [BFGS].
- Parameters
h0 (
float
) – Initial guess of kernel bandwidth in pixels.method (
str
) –reference_image update algorithm. The following keyword values are allowed:
KerReg : Kernel regression algorithm.
LOWESS : Local weighted linear regression.
beta (
float
) – Exponential decay coefficient for inverse Hessian matrix.epsilon (
float
) – The step size used in the estimation of the fitness gradient.maxiter (
int
) – Maximum number of iterations in the minimization loop.gtol (
float
) – Gradient norm must be less than gtol before successful termination.verbose (
bool
) – Print convergence message if True.
- Return type
- Returns
Optimal kernel bandwidth in pixels.
References
- BFGS
S. Wright, J. Nocedal et al., “Numerical optimization,” Springer Sci. 35, 7 (1999).
- ref_indices()#
Return an array of reference profile pixel indices.
- Return type
- Returns
An array of reference profile pixel indices.
- test_train_split(test_ratio=0.1)#
Update test / train subsets split. Required to calculate the Cross-validation error metric.
- Parameters
test_ratio (
float
) – Ratio between the size of the test subset and the whole dataset.- Return type
- Returns
A new
SpeckleTracking
object with a new test / train subsets split.
See also
SpeckleTracking.CV()
: Full details on the Cross-validation error metric.
- train(search_window, h0, blur=0.0, n_iter=30, f_tol=1e-08, ref_method='KerReg', pm_method='rsearch', pm_args={}, options={}, verbose=True)#
Perform iterative Robust Speckle Tracking update. The reconstruction cycle consists of: (i) generating the reference image (see
SpeckleTracking.update_reference()
); (ii) updating the pixel mapping between a stack of frames and the generated reference image (seeSpeckleTracking.update_pixel_map()
); (iii) updating the sample translation vectors (if needed, seeSpeckleTracking.update_translations()
); and (iv) calculating the mean Huber error (seeSpeckleTracking.update_error()
). The kernel bandwidth in the reference update is kept fixed during the iterative update procedure.- Parameters
search_window (
Tuple
[float
,float
,float
]) –A tuple of three elements (‘sw_y’, ‘sw_x’, ‘sw_s’). The elements are the following:
sw_y : Search window size in pixels along the horizontal detector axis.
sw_x : Search window size in pixels along the vertical detector axis.
sw_s : Search window size of the Huber scaling map. Given as a ratio (0.0 - 1.0) relative to the scaling map value before the update.
h0 (
float
) – Smoothing kernel bandwidth used in reference_image regression. The value is given in pixels.blur (
float
) – Smoothing kernel bandwidth used in pixel_map regularisation. The value is given in pixels.n_iter (
int
) – Maximum number of iterations.f_tol (
float
) – Tolerance for termination by the change of the average error. The iteration stops when(f^k - f^{k+1})/max{|f^k|,|f^{k+1}|,1} <= ftol
.ref_method (
str
) –reference_image update algorithm. The following keyword values are allowed:
KerReg : Kernel regression algorithm.
LOWESS : Local weighted linear regression.
pm_method (
str
) –pixel_map update algorithm. The following keyword values are allowed:
gsearch : Grid search algorithm.
rsearch : Random search algorithm.
de : Differential evolution algorithm.
pm_args (
Dict
[str
,Union
[bool
,int
,float
,str
]]) –Pixel map update options. Accepts the following keyword arguments:
integrate : Ensure that the updated pixel map is irrotational by integrating and taking the derivative. False by default.
grid_size : Grid size along one of the detector axes for the ‘gsearch’ method. The grid shape is then (grid_size, grid_size). The default value is
int(0.5 * (sw_x + sw_y))
.n_trials : Number of points generated at each pixel in the detector grid for the ‘rsearch’ method. The default value is
int(0.25 * (sw_x + sw_y)**2)
.n_iter : The maximum number of generations over which the entire population is evolved in the ‘de’ method. The default value is 5.
pop_size : The total population size in the ‘de’ method. Must be greater or equal to 4. The default value is
max(int(0.25 * (sw_x + sw_y)**2) / n_iter, 4)
.mutation : The mutation constant in the ‘de’ method. It should be in the range [0, 2]. The default value is 0.75.
recombination : The recombination constant in the ‘de’ method, should be in the range [0, 1]. The default value is 0.7.
seed : Specify seed for the random number generation. Generated automatically if not provided.
options (
Dict
[str
,Union
[bool
,float
,str
]]) –Extra options. Accepts the following keyword arguments:
momentum : Momentum used in the next error calculation. Helps to smooth out the change of error.
update_translations : Update sample pixel translations if True. The default value is False.
return_extra : Return errors at each iteration if True. The default value is False.
verbose (
bool
) – Set verbosity of the computation process.
- Return type
- Returns
A tuple of two items (st_obj, errors). The elements of the tuple are as follows:
st_obj : A new
SpeckleTracking
object with the updated pixel_map and reference_image. di_pix and dj_pix are also updated if ‘update_translations’ in options is True.errors : List of average error values for each iteration. Only if ‘return_extra’ in options is True.
- train_adapt(search_window, h0, blur=0.0, n_iter=30, f_tol=1e-08, ref_method='KerReg', pm_method='rsearch', pm_args={}, options={}, verbose=True)#
Perform adaptive iterative Robust Speckle Tracking update. The reconstruction cycle consists of: (i) estimating an optimal kernel bandwidth for the reference image estimate (see
SpeckleTracking.find_hopt()
); (ii) generating the reference image (seeSpeckleTracking.update_reference()
); (iii) updating the pixel mapping between a stack of frames and the generated reference image (seeSpeckleTracking.update_pixel_map()
); (iv) updating the sample translation vectors (if needed, seeSpeckleTracking.update_translations()
); and (v) calculating the mean Huber error (seeSpeckleTracking.update_errors()
).- Parameters
search_window (
Tuple
[float
,float
,float
]) –A tuple of three elements (‘sw_y’, ‘sw_x’, ‘sw_s’). The elements are the following:
sw_y : Search window size in pixels along the horizontal detector axis.
sw_x : Search window size in pixels along the vertical detector axis.
sw_s : Search window size of the Huber scaling map. Given as a ratio (0.0 - 1.0) relative to the scaling map value before the update.
h0 (
float
) – Smoothing kernel bandwidth used in reference_image estimation. The value is given in pixels.blur (
float
) – Smoothing kernel bandwidth used in pixel_map regularisation. The value is given in pixels.n_iter (
int
) – Maximum number of iterations.f_tol (
float
) – Tolerance for termination by the change of the average error. The iteration stops when(error^k - error^{k+1})/max{|error^k|, |error^{k+1}|, 1} <= f_tol
.ref_method (
str
) –reference_image update algorithm. The following keyword values are allowed:
KerReg : Kernel regression algorithm.
LOWESS : Local weighted linear regression.
pm_method (
str
) –pixel_map update algorithm. The following keyword values are allowed:
gsearch : Grid search algorithm.
rsearch : Random search algorithm.
de : Differential evolution algorithm.
pm_args (
Dict
[str
,Union
[bool
,int
,float
,str
]]) –Pixel map update options. Accepts the following keyword arguments:
integrate : Ensure that the updated pixel map is irrotational by integrating and taking the derivative. False by default.
grid_size : Grid size along one of the detector axes for the ‘gsearch’ method. The grid shape is then (grid_size, grid_size). The default value is
int(0.5 * (sw_x + sw_y))
.n_trials : Number of points generated at each pixel in the detector grid for the ‘rsearch’ method. The default value is
int(0.25 * (sw_x + sw_y)**2)
.n_iter : The maximum number of generations over which the entire population is evolved in the ‘de’ method. The default value is 5.
pop_size : The total population size in the ‘de’ method. Must be greater or equal to 4. The default value is
max(int(0.25 * (sw_x + sw_y)**2) / n_iter, 4)
.mutation : The mutation constant in the ‘de’ method. It should be in the range [0, 2]. The default value is 0.75.
recombination : The recombination constant in the ‘de’ method, should be in the range [0, 1]. The default value is 0.7.
seed : Specify seed for the random number generation. Generated automatically if not provided.
options (
Dict
[str
,Union
[bool
,float
,str
]]) –Extra options. Accepts the following keyword arguments:
beta : Exponential decay coefficient for inverse Hessian matrix. The default value is 0.9.
epsilon : Increment to h0 to use for determining the function gradient for h0 update algorithm. The default value is 1.4901161193847656e-08.
maxiter : Maximum number of iterations in the line search at the optimal kernel bandwidth update. The default value is 10.
momentum : Momentum used in the next error calculation. Helps to smooth out the change of error. The default value is 0.0.
update_translations : Update sample pixel translations if True. The default value is False.
return_extra : Return errors and h0 array if True. The default value is False.
verbose (
bool
) – Set verbosity of the computation process.
- Return type
- Returns
A tuple of two items (‘st_obj’, ‘extra’). The elements of the tuple are as follows:
st_obj : A new
SpeckleTracking
object with the updated pixel_map and reference_image. di_pix and dj_pix are also updated if update_translations is True.extra: A dictionary with the given parameters:
errors : List of average error values for each iteration.
hvals : List of kernel bandwidths for each iteration.
Only if return_extra is True.
- update_errors()#
Return a new
SpeckleTracking
object with the updated mean Huber error metric error.- Return type
- Returns
A new
SpeckleTracking
object with the updated error.- Raises
AttributeError – If reference_image was not generated beforehand.
See also
pyrost.bin.ref_errors()
: Full details of the reference update error metric.pyrost.bin.pm_errors()
: Full details of the pixel mapping update error metric.
- update_pixel_map(search_window, blur=0.0, integrate=False, method='gsearch', extra_args={})#
Return a new
SpeckleTracking
object with the updated pixel mapping (pixel_map) and Huber scale mapping (scale_map). The update is performed with the adaptive Huber regression [HUBER]. The Huber loss function is minimized by the non-gradient methods enlisted in the method argument.- Parameters
search_window (
Tuple
[float
,float
,float
]) –A tuple of three elements (‘sw_y’, ‘sw_x’, ‘sw_s’). The elements are the following:
sw_y : Search window size in pixels along the horizontal detector axis.
sw_x : Search window size in pixels along the vertical detector axis.
sw_s : Search window size of the Huber scaling map. Given as a ratio (0.0 - 1.0) relative to the scaling map value before the update.
blur (
float
) – Smoothing kernel bandwidth used in pixel_map regularisation. The value is given in pixels.integrate (
bool
) – Ensure that the updated pixel map is irrotational by integrating and taking the derivative.method (
str
) –pixel_map update algorithm. The following keyword values are allowed:
gsearch : Grid search algorithm.
rsearch : Random search algorithm.
de : Differential evolution algorithm.
extra_args (
Dict
[str
,Union
[int
,float
]]) –Extra arguments for pixel map update methods. Accepts the following keyword arguments:
grid_size : Grid size along one of the detector axes for the ‘gsearch’ method. The grid shape is then (grid_size, grid_size). The default value is
int(0.5 * (sw_x + sw_y))
.n_trials : Number of points generated at each pixel in the detector grid for the ‘rsearch’ method. The default value is
int(0.25 * (sw_x + sw_y)**2)
.n_iter : The maximum number of generations over which the entire population is evolved in the ‘de’ method. The default value is 5.
pop_size : The total population size in the ‘de’ method. Must be greater or equal to 4. The default value is
max(int(0.25 * (sw_x + sw_y)**2) / n_iter, 4)
.mutation : The mutation constant in the ‘de’ method. It should be in the range [0, 2]. The default value is 0.75.
recombination : The recombination constant in the ‘de’ method, should be in the range [0, 1]. The default value is 0.7.
seed : Specify seed for the random number generation. Generated automatically if not provided.
- Raises
AttributeError – If reference_image was not generated beforehand.
ValueError – If method keyword value is not valid.
- Return type
- Returns
A new
SpeckleTracking
object with the updated pixel_map and scale_map.
Notes
The pixel mapping update is carried out separately at each pixel \({i, j}\) in the detector grid by minimizing the Huber error metric given by:
\[\varepsilon_{ij}(\delta i, \delta j, s) = \frac{1}{N} \sum_{n = 1}^N \left[ s + \mathcal{H}_{1.35} \left( \frac{I_{nij} - W_{ij} I_\text{ref}(f_x i - u^x_{ij} - \delta i + \Delta i_n, f_y i - u^y_{ij} - \delta j + \Delta j_n)}{s} \right) s \right],\]where \(I_\text{ref}\) is the reference profile, \(u^x_{ij}, u^y_{ij}\) are the horizontal and vertical components of the pixel mapping, \(\Delta i_n, \Delta j_n\) are the sample translation along the horizontal and vertical axes in pixels, \(I_{nij}\) are the measured stack of frames, and W_{ij} is the white-field.
See also
pyrost.bin.pm_gsearch()
: Full details of the grid search update method.pyrost.bin.pm_rsearch()
: Full details of the random search update method.pyrost.bin.pm_devolution()
: Full details of the differential evolution update method.
References
- update_reference(hval, method='KerReg')#
Return a new
SpeckleTracking
object with the updated reference_image. The reference profile is estimated either by Kernel regression (‘KerReg’) [KerReg] or Local weighted linearegressin (‘LOWESS’) [LOWESS].- Parameters
- Raises
ValueError – If method keyword value is not valid.
- Return type
- Returns
A new
SpeckleTracking
object with the updated reference_image.
Notes
The reference profile estimator is given by:
\[I_\text{ref}(f_x i, f_y j) = \frac{\sum_n \sum_{i^{\prime}} \sum_{j^{\prime}} K(f_x i - u^x_{i^{\prime}j^{\prime}} + \Delta i_n, f_y j - u^y_{i^{\prime}j^{\prime}} + \Delta j_n, h) \; W_{i^{\prime}j^{\prime}} I_{n i^{\prime} j^{\prime}}} {\sum_n \sum_{i^{\prime}} \sum_{j^{\prime}} K(f_x i - u^x_{i^{\prime}j^{\prime}} + \Delta i_n, f_y j - u^y_{i^{\prime}j^{\prime}} + \Delta j_n, h) \; W^2_{i^{\prime} j^{\prime}}},\]where \(K(i, j, h) = \exp(\frac{i\:\delta u + j\:\delta v}{h}) / \sqrt{2 \pi}\) is the Gaussian kernel, \(u^x_{ij}, u^y_{ij}\) are the horizontal and vertical components of the pixel mapping, \(\Delta i_n, \Delta j_n\) are the sample translation along the horizontal and vertical axes in pixels, \(I_{nij}\) are the measured stack of frames, and W_{ij} is the white-field.
See also
pyrost.bin.make_reference()
Full details of the reference_imageupdate algorithm.
References
- update_translations(sw_x, sw_y=0.0, blur=0.0)#
Return a new
SpeckleTracking
object with the updated sample pixel translations (di_pix, dj_pix). The update is performed with the adaptive Huber regression [HUBER]. The Huber loss function is minimized by the grid search algorithm.- Parameters
- Return type
- Returns
A new
SpeckleTracking
object with the updated di_pix, dj_pix.- Raises
AttributeError – If reference_image was not generated beforehand.
See also
pyrost.bin.tr_gsearch()
Full details of the sample translationsupdate algorithm.