STData#
- class pyrost.STData(input_file, transform=None, basis_vectors=None, data=None, defocus_x=None, defocus_y=None, distance=None, frames=None, output_file=None, phase=None, pixel_aberrations=None, reference_image=None, scale_map=None, translations=None, wavelength=None, whitefields=None, x_pixel_size=None, y_pixel_size=None, good_frames=None, mask=None, num_threads=2, pixel_translations=None, whitefield=None)#
Speckle tracking data container class. Needs a
pyrost.CXIStore
file handler. Provides an interface to work with the data and contains a suite of tools for the R-PXST data processing pipeline. Also provides an interface to load from a file and save to a file any of the data attributes. The data frames can be tranformed using any of thepyrost.Transform
classes.- Parameters
input_file (
CXIStore
) – HDF5 or CXI file handler of input files.kwargs – Dictionary of the necessary and optional data attributes specified in
pyrost.STData
notes. All the necessary attributes must be provided
- Raises
ValueError – If any of the necessary attributes specified in
pyrost.STData
notes have not been provided.
Notes
Necessary attributes:
basis_vectors : Detector basis vectors
data : Measured intensity frames.
distance : Sample-to-detector distance [m].
frames : List of frame indices.
translations : Sample’s translations [m].
wavelength : Incoming beam’s wavelength [m].
x_pixel_size : Pixel’s size along the horizontal detector axis [m].
y_pixel_size : Pixel’s size along the vertical detector axis [m].
Optional attributes:
defocus_x : Defocus distance for the horizontal detector axis [m].
defocus_y : Defocus distance for the vertical detector axis [m].
good_frames : An array of good frames’ indices.
mask : Bad pixels mask.
num_threads : Number of threads used in computations.
output_file : Output file handler.
phase : Phase profile of lens’ aberrations.
pixel_aberrations : Lens’ aberrations along the horizontal and vertical axes in pixels.
pixel_translations : Sample’s translations in the detector’s plane in pixels.
reference_image : The unabberated reference image of the sample.
scale_map : Huber scale map.
whitefield : Measured frames’ white-field.
whitefields : Set of dynamic white-fields for each of the measured images.
- clear(attributes=None)#
Clear the data inside the container.
- contents()#
Return a list of the attributes stored in the container that are initialised.
- defocus_sweep(defoci_x, defoci_y=None, size=51, hval=None, extra_args={}, return_extra=False, verbose=True)#
Calculate a set of reference images for each defocus in defoci and return an average R-characteristic of an image (the higher the value the sharper reference image is). The kernel bandwidth hval is automatically estimated by default. Return the intermediate results if return_extra is True.
- Parameters
defoci_x (
ndarray
) – Array of defocus distances along the horizontal detector axis [m].defoci_y (
Optional
[ndarray
]) – Array of defocus distances along the vertical detector axis [m].hval (
Optional
[float
]) – Kernel bandwidth in pixels for the reference image update. Estimated withpyrost.SpeckleTracking.find_hopt()
for an average defocus value if None.size (
int
) – Local variance filter size in pixels.extra_args (
Dict
[str
,Union
[bool
,float
,str
]]) –Extra arguments parser to the
STData.get_st()
andSpeckleTracking.update_reference()
methods. The following keyword values are allowed:ds_y : Reference image sampling interval in pixels along the horizontal axis. The default value is 1.0.
ds_x : Reference image sampling interval in pixels along the vertical axis. The default value is 1.0.
aberrations : Add pixel_aberrations to pixel_map of
SpeckleTracking
object if it’s True. The default value is False.ff_correction : Apply dynamic flatfield correction if it’s True. The default value is False.
ref_method : Choose the reference image update algorithm. The following keyword values are allowed:
KerReg : Kernel regression algorithm.
LOWESS : Local weighted linear regression.
The default value is ‘KerReg’.
return_extra (
bool
) – Return a dictionary with the intermediate results if True.verbose (
bool
) – Set the verbosity of the process.
- Return type
- Returns
A tuple of two items (‘r_vals’, ‘extra’). The elements are as follows:
r_vals : Array of the average values of reference_image gradients squared.
extra : Dictionary with the intermediate results. Only if return_extra is True. Contains the following data:
reference_image : The generated set of reference profiles.
r_images : The set of local variance images of reference profiles.
Notes
R-characteristic is called a local variance and is given by:
\[R[i, j] = \frac{\sum_{i^{\prime} = -N / 2}^{N / 2} \sum_{j^{\prime} = -N / 2}^{N / 2} (I[i - i^{\prime}, j - j^{\prime}] - \bar{I}[i, j])^2}{\bar{I}^2[i, j]},\]where \(\bar{I}[i, j]\) is a local mean and defined as follows:
\[\bar{I}[i, j] = \frac{1}{N^2} \sum_{i^{\prime} = -N / 2}^{N / 2} \sum_{j^{\prime} = -N / 2}^{N / 2} I[i - i^{\prime}, j - j^{\prime}]\]See also
pyrost.SpeckleTracking.update_reference()
: reference image update algorithm.
- fit_phase(center=0, axis=1, max_order=2, xtol=1e-14, ftol=1e-14, loss='cauchy')#
Fit pixel_aberrations with the polynomial function using nonlinear least-squares algorithm. The function uses least-squares algorithm from
scipy.optimize.least_squares()
.- Parameters
center (
int
) – Index of the zerro scattering angle or direct beam pixel.axis (
int
) – Axis along which pixel_aberrations is fitted.max_order (
int
) – Maximum order of the polynomial model function.xtol (
float
) – Tolerance for termination by the change of the independent variables.ftol (
float
) – Tolerance for termination by the change of the cost function.loss (
str
) –Determines the loss function. The following keyword values are allowed:
linear :
rho(z) = z
. Gives a standard least-squares problem.soft_l1 :
rho(z) = 2 * ((1 + z)**0.5 - 1)
. The smooth approximation of l1 (absolute value) loss. Usually a good choice for robust least squares.huber :
rho(z) = z if z <= 1 else 2*z**0.5 - 1
. Works similarly to ‘soft_l1’.cauchy (default) :
rho(z) = ln(1 + z)
. Severely weakens outliers influence, but may cause difficulties in optimization process.arctan :
rho(z) = arctan(z)
. Limits a maximum loss on a single residual, has properties similar to ‘cauchy’.
- Return type
- Returns
A dictionary with the model fit information. The following fields are contained:
c_3 : Third order aberrations coefficient [rad / mrad^3].
c_4 : Fourth order aberrations coefficient [rad / mrad^4].
fit : Array of the polynomial function coefficients of the pixel aberrations fit.
ph_fit : Array of the polynomial function coefficients of the phase aberrations fit.
rel_err : Vector of relative errors of the fit coefficients.
r_sq :
R**2
goodness of fit.
See also
pyrost.AberrationsFit.fit()
: Full details of the aberrations fitting algorithm.
- get(attr, value=None)#
Retrieve a dataset, return
value
if the attribute is not found.
- get_fit(center=0, axis=1)#
Return an
AberrationsFit
object for parametric regression of the lens’ aberrations profile. Raises an error if ‘defocus_x’ or ‘defocus_y’ is not defined.- Parameters
- Raises
ValueError – If ‘defocus_x’ or ‘defocus_y’ is not defined in the container.
- Return type
- Returns
An instance of
AberrationsFit
class.
- get_pca()#
Perform the Principal Component Analysis [PCA] of the measured data and return a set of eigen flatfields (EFF).
- Return type
- Returns
A tuple of (‘cor_data’, ‘effs’, ‘eig_vals’). The elements are as follows:
cor_data : Background corrected stack of measured frames.
effs : Set of eigen flat-fields.
eig_vals : Corresponding eigen values for each of the eigen flat-fields.
References
- get_st(ds_y=1.0, ds_x=1.0, test_ratio=0.1, aberrations=False, ff_correction=False)#
Return
SpeckleTracking
object derived from the container. Return None if defocus_x or defocus_y doesn’t exist in the container.- Parameters
ds_y (
float
) – Reference image sampling interval in pixels along the vertical axis.ds_x (
float
) – Reference image sampling interval in pixels along the horizontal axis.test_ratio (
float
) – Ratio between the size of the test subset and the whole dataset.aberrations (
bool
) – Add pixel_aberrations to pixel_map ofSpeckleTracking
object if it’s True.ff_correction (
bool
) – Apply dynamic flatfield correction if it’s True.
- Return type
- Returns
An instance of
SpeckleTracking
derived from the container. None if defocus_x or defocus_y are not defined.
- import_st(st_obj)#
Update pixel_aberrations, phase, reference_image, and scale_map based on the data from st_obj object. st_obj must be derived from this data container, an error is raised otherwise.
- Parameters
st_obj (
SpeckleTracking
) –SpeckleTracking
object derived from this data container.- Raises
ValueError – If st_obj wasn’t derived from this data container.
- integrate_data(axis=0)#
Return a new
STData
object with the data summed over the axis. Clear all the 2D and 3D data attributes inside the container.
- items()#
Return (key, value) pairs of the datasets stored in the container.
- Return type
- Returns
(key, value) pairs of the datasets stored in the container.
- keys()#
Return a list of the attributes available in the container.
- load(attributes=None, idxs=None, processes=1, verbose=True)#
Load data attributes from the input files in files file handler object.
- Parameters
attributes (
Union
[str
,List
[str
],None
]) – List of attributes to load. Loads all the data attributes contained in the file(s) by default.idxs (
Union
[int
,slice
,ndarray
,List
[int
],Tuple
[int
],None
]) – List of frame indices to load.processes (
int
) – Number of parallel workers used during the loading.verbose (
bool
) – Set the verbosity of the loading process.
- Raises
ValueError – If attribute is not existing in the input file(s).
ValueError – If attribute is invalid.
- Return type
- Returns
New
STData
object with the attributes loaded.
- mask_frames(frames=None)#
Return a new
STData
object with the updated good frames mask. Mask empty frames by default.
- pixel_map(dtype=<class 'float64'>)#
Return a preliminary pixel mapping.
- replace(**kwargs)#
Return a new
pyrost.STData
container with replaced data.- Parameters
kwargs (Any) – Replaced attributes.
- Return type
- Returns
A new
pyrost.STData
container.
- save(attributes=None, apply_transform=False, mode='append', idxs=None)#
Save data arrays of the data attributes contained in the container to an output file.
- Parameters
attributes (
Union
[str
,List
[str
],None
]) – List of attributes to save. Saves all the data attributes contained in the container by default.apply_transform (
bool
) – Apply transform to the data arrays if True.mode (
str
) –Writing modes. The following keyword values are allowed:
append : Append the data array to already existing dataset.
insert : Insert the data under the given indices idxs.
overwrite : Overwrite the existing dataset.
idxs (
Union
[int
,slice
,ndarray
,List
[int
],Tuple
[int
],None
]) – Indices where the data is saved. Used only ifmode
is set to ‘insert’.
- Raises
ValueError – If the
output_file
is not defined inside the container.
- update_defocus(defocus_x, defocus_y=None)#
Return a new
STData
object with the updated defocus distances defocus_x and defocus_y for the horizontal and vertical detector axes accordingly. Update pixel_translations based on the new defocus distances.- Parameters
- Return type
- Returns
New
STData
object with the updated defocus_y, defocus_x, and pixel_translations.
- update_mask(method='perc-bad', pmin=0.0, pmax=99.99, vmin=0, vmax=65535, update='reset')#
Return a new
STData
object with the updated bad pixels mask.- Parameters
method (
str
) –Bad pixels masking methods. The following keyword values are allowed:
’no-bad’ (default) : No bad pixels.
’range-bad’ : Mask the pixels which values lie outside of (vmin, vmax) range.
’perc-bad’ : Mask the pixels which values lie outside of the (pmin, pmax) percentiles.
vmin (
int
) – Lower intensity bound of ‘range-bad’ masking method.vmax (
int
) – Upper intensity bound of ‘range-bad’ masking method.pmin (
float
) – Lower percentage bound of ‘perc-bad’ masking method.pmax (
float
) – Upper percentage bound of ‘perc-bad’ masking method.update (
str
) – Multiply the new mask and the old one if ‘multiply’, use the new one if ‘reset’.
- Raises
ValueError – If there is no
data
inside the container.ValueError – If
method
keyword is invalid.ValueError – If
update
keyword is invalid.ValueError – If
vmin
is larger thanvmax
.ValueError – If
pmin
is larger thanpmax
.
- Return type
- Returns
New
STData
object with the updatedmask
.
- update_whitefields(method='median', size=11, cor_data=None, effs=None)#
Return a new
STData
object with a new set of dynamic whitefields. A set of whitefields are generated by the dint of median filtering or Principal Component Analysis [PCA].- Parameters
method (
str
) –Method to generate a set of dynamic white-fields. The following keyword values are allowed:
median : Median data along the first axis.
pca : Generate a set of dynamic white-fields based on eigen flatfields effs from the PCA. effs can be obtained with
STData.get_pca()
method.
size (
int
) – Size of the filter window in pixels used for the ‘median’ generation method.cor_data (
Optional
[ndarray
]) – Background corrected stack of measured frames.effs (
Optional
[ndarray
]) – Set of Eigen flatfields used for the ‘pca’ generation method.
- Raises
ValueError – If the method keyword is invalid.
AttributeError – If the whitefield is absent in the
STData
container when using the ‘pca’ generation method.ValueError – If effs were not provided when using the ‘pca’ generation method.
- Return type
- Returns
New
STData
object with the updated whitefields.
See also
pyrost.STData.get_pca()
: Method to generate eigen flatfields.
- values()#
Return the attributes’ data stored in the container.
- Return type
- Returns
List of data stored in the container.