Usage¶
Dataset format¶
Adorym reads raw data contained an HDF5 file. The diffraction images
should be
stored in the exchange/data dataset as a 4D array, with a shape of
[n_rotation_angles, n_diffraction_spots, image_size_y, image_size_x].
In a large part, Adorym is blind to the type of experiment, which
means
there no need to explicitly tell it the imaging technique used to
generate
the dataset. For imaging data collected from only one angle,
n_rotation_angles = 1.
For full-field imaging without scanning, n_diffraction_spots = 1.
For
2D imaging, set the last dimension of the object size to 1 (this will
be
introduced further below).
Experimental metadata including beam energy, probe position, and pixel
size, may also be stored in the HDF5, but they can also be provided
individually
as arguments to the function reconstruct_ptychography. When these
arguments
are provided, Adorym uses the arguments rather than reads the metadata
from
the HDF5.
The following is the full structure of the HDf5:
data.h5
|___ exchange
| |___ data: float, 4D array
| [n_rotation_angles, n_diffraction_spots, image_size_y, image_size_x]
|
|___ metadata
|___ energy_ev: scalar, float. Beam energy in eV
|___ probe_pos_px: float, [n_diffraction_spots, 2].
| Probe positions (y, x) in pixel.
|___ psize_cm: scalar, float. Sample-plane pixel size in cm.
|___ free_prop_cm: (optional) scalar or array
| Distance between sample exiting plane and detector.
| For far-field propagation, do not include this item.
|___ slice_pos_cm: (optional) float, 1D array
Position of each slice in sparse multislice ptychography. Starts from 0.
Parameter settings in main function¶
The scripts in demos and tests supply the
reconstruct_ptychography
with parameters listed as a Python dictionary. You may find the
docstrings
of the function helpful, but here lists a collection of the most
crucial
parameters:
Backend¶
| Arg name | Type | Default | Description |
|---|---|---|---|
backend |
String | autograd |
Select 'pytorch' or 'autograd'. Both can be used as the automatic differentiation engine, but only the PyTorch backend supports GPU computation. Some features are only supported through PyTorch, including affine transformation refinement and object tilt-angle refinement. |
Raw data and experimental parameters¶
| Arg name | Type | Default | Description |
|---|---|---|---|
fname |
String | Name of the HDF5 containing raw data. Put only the basename here; any path predix should go to save_path. Some features are only supported through PyTorch, including affine transformation refinement and object tilt-angle refinement. |
|
obj_size |
Array of Int | [L_y, L_x, L_z]. The size of the object function (i.e., the unknowns) in pixels. L_y is the size in the vertical direction, while L_x and L_z refer to sizes on the horizontal plane. For 2D reconstruction, set L_z to 1. For 3D reconstruction, it is strongly recommended to keep L_x == L_z. For doing sparsely spaced multislice tomography (i.e., when the number of slices along beam axis is much less than the number of lateral pixels), the best practice is to set binning to a larger value, instead of using a small L_z. |
|
probe_pos |
Array of Float | None |
[n_diffraction_spots, 2]. Probe positions in a scanning-type experiment in pixel in the object frame (i.e., real-unit probe positions divided by sample plane pixel size). Default to None. If None, the program will attempt to get the value from HDF5. The positions will be interpreted as the top-left corner of the probe array in object frame. For single-shot experiments, set probe_pos to [[0, 0]]. |
theta_st |
Float | 0 | Starting rotation angle in radian. Default to 0. |
theta_end |
Float | PI |
End rotation angle in radian. For single angle data, set this the same as theta_st. |
n_theta |
Int | None |
Number of rotation angles. If None, the number will be inferred from the shape of the HDF5 dataset. |
theta_downsample |
Int | None |
By how many times the raw data should be downsampled in rotation angles. |
energy_ev |
Float | None |
X-ray beam energy in eV. If None, the program will attempt to get the value from HDF5. |
psize_cm |
Float | None |
Lateral pixel size at sample plane in cm. If None, the program will attempt to get the value from HDF5. If axial pixel size is different, use slice_pos_cm. |
free_prop_cm |
Float | None |
The distance between sample and detector in cm. For far-field imaging, set it to None or 'inf', so that the programs uses Fraunhofer approximation. For near-field imaging, this value is assumed to be the propagation distance in a plane-wave illuminated experiment; if the illumination is a spherical wave generated by a point source, use the effective distance given by Fresnel scaling theorem: ``z_eff = z1 * z2 / (z1 + z2)``. |
raw_data_type |
String | 'intensity' |
Choose from 'intensity' or 'magnitude'. This informs the optimizer the type of raw data contained in the HDF5, and determines whether the measured data should be square-rooted when calculating loss. For conventional tomography with ``pure_propjection=True`` and ``is_minus_logged=True``, this must be ``magnitude``! |
is_minus_logged |
Boolean | False |
Whether the raw projection data have been minus-logged. This is usually used in conventional tomography. If True, forward model will return a simple summation of beta along the beam axis. |
slice_pos_cm |
Array of Float | None |
Position of each slice in sparse multislice ptychography. Starts from 0. If None, the program will attempt to get the value from HDF5. |
Reconstruction parameters¶
| Arg name | Type | Default | Description |
|---|---|---|---|
n_epochs |
Int | 'auto' |
Number of epochs to run. An epoch refers to a cycle during which all diffraction data are processed. Set it to 'auto' to automatically stops the reconstruction when the reduction rate of loss falls below crit_conv_rate. This option is not recommended especially for noisy data due to the possibility of fake positives. The best practice so far is to set n_epochs to a sufficiently large value and observe the loss curve and reconstruction output until satisfactory results are obtained. |
crit_conv_rate |
Float | 0.03 | If the reduction rate of loss at the current epoch in regards to the previous one is below this value, convergence is assumed to be reached and the reconstruction process stops. |
max_epochs |
Int | 200 | When n_epochs is set to 'auto', the program will stop regardless of the loss reduction rate once this number of epochs have been run. |
regularizers |
List | None |
A list of Regularizer objects. Alternatively, you can specify regularizers by keeping this argument as None and supply alpha_d, alpha_b, and gamma instead. |
alpha_d |
Float | 0 | Weight applied to l1-norm of the delta (or real) part of the object function, depending on the setting of unknown_type. The full loss function is in the form of L = D(f(x), y0) + alpha_d * |x_d|_1 + alpha_b * |x_b|_1 + gamma * TV(x). Ignored when regularizers is not None. |
alpha_b |
Float | 0 | Weight applied to l1-norm of the beta (or imaginary) part of the object function. Ignored when regularizers is not None. |
gamma |
Float | 0 | Weight applied to total variation of the object function. Ignored when regularizers is not None. |
minibatch_size |
Int | 1 | The number of diffraction spots to be processed at a time. When multi-processing, this is the number of diffraction spots processed by each rank. |
multiscale_level |
Int | 1 | Number of levels for multi-scale progressive reconstruction. This feature is still experimental. |
n_epoch_final_pass |
Int | None | If multiscale_level is larger than 1, this parameter sets the number of epochs for the last (full-resolution) pass. |
initial_guess |
List of Arrays | None | The initial guess of the object function in the form of [obj_delta, obj_beta] when unknown_type is delta_beta, or [obj_mag, obj_phase] when unknown_type is real_imag. The arrays must have the same size as specified by obj_size. |
random_guess_means_sigmas |
List of Floats | (8.7e-7, 5.1e-8, 1e-7, 1e-8) |
When initial_guess is None, the object function will be initialized usin Gaussian randoms. This argument provides the Gaussian parameters in the format of (mean_delta, mean_beta, sigma_delta, sigma_beta) or (mean_mag, mean_phase, sigma_mag, sigma_phase), depending on the setting of unknwon_type. |
n_batch_per_update |
Int | 1 | The number of minibatches to accumulate before the object is updated. Ignored when update_scheme is per angle. |
reweighted_l1 |
Bool | False |
If True and alpha_d != 0, the program uses reweighted l1-norm to regularize the object (see Candès, E. J., Wakin, M. B. & Boyd, S. P. Enhancing Sparsity by Reweighted ℓ1 Minimization. Journal of Fourier Analysis and Applications 14, (2008). ) |
interpolation |
String | 'bilinear' |
Interpolation method for rotation. |
update_scheme |
String | 'immediate' |
Choose from 'immediate' or 'per angle'. If 'immediate', the object function is updated immedaitely after each minibatch is done. If 'per angle', updated is performed only after all diffraction patterns from the current rotation angle are processed. If shared_file_object is on, the 'per angle' mode is used regardless of this setting. |
unknown_type |
String | 'delta_beta' |
Choose from delta_beta and real_imag. If set to delta_beta, the program treats the unknowns as the delta and beta parts in the complex refractive indices of the object, n = 1-delta-i*beta. In this case, modulation to the wavefield by each slice of the object will be done as wavefield * exp(-i*k*n*z). If set to real_imag, the unknowns are treated as the real and imaginary part of a multiplicative object function, where the modulation is done as wavefield * (obj_real + i * obj_imag). Using delta_beta can help overcome mild phase wrapping, while using real_imag generally leads to better numerical robustness. |
randomize_probe_pos |
Bool | False | Whether to randomize diffraction spots on each viewing angle when there are more than 1 of them. Recommended to be True for 2D ptychography. |
common_probe_pos |
Bool | True | Whether the number and position of tiles are the same for all viewing angles. If False, the tile positions for each angle should be provided in the HDF5 as ‘metadata/probe_pos_px_’. The main dataset remains as a 4D array, where the size of the second axis is determined by the angle that has the most tiles. |
beamstop |
Tensor | None | A 2D boolean mask, if not None. When supplied, only unmasked detector pixels are used for computing the mismatch term. |
Object optimizer options¶
| Arg name | Type | Default | Description |
|---|---|---|---|
optimize_object |
Bool | True |
Keep True in most cases. Setting to False forbids the object from being updated using gradients, which might be desirable when you just want to refine parameters for other reconstruction algorithms. |
optimizer |
adorym.Optimizer or String |
'adam' |
Either a predeclared adorym.Optimizer class, or choose from 'adam', 'gd' (steepest gradient descent), 'momentum', or 'cg'. You may also try 'curveball' but it is still experimental and supports only data parallelism mode. |
learning_rate |
Float | 1e-5 |
Learning rate, or step size of the chosen optimizer for the object function. Ignored if optimizer is 'curveball'. |
optimizer_batch_number_increment |
String | 'angle' |
Applies to optimizers that use the current batch number for calculation, such as Adam. If 'angle', batch number passed to optimizer increments after each angle. This is recommended for 2D reconstruction. If 'batch', it increases after each batch. This is recommended for 3D reconstruction. If distribution_mode is not None, 'batch' behaves the same as 'angle'. |
Finite support constraint¶
| Arg name | Type | Default | Description |
|---|---|---|---|
finite_support_mask_path |
String | None |
The path to the TIFF file storing the finite support mask. In general, this is needed only for single-shot CDI and holography. |
shrink_cycle |
Int | None |
For every how many minibatches should the finite support mask be shrink-wrapped. Use None to disable shrink-wrap. Useful only when finite_support_mask_path is not None. |
'shrink_threshold' |
Float | 1e-9 |
Threshold for shrink-wrapping. Useful only when finite_support_mask_path is not None. |
Object contraints¶
| Arg name | Type | Default | Description |
|---|---|---|---|
object_type |
String | 'normal' |
Choose from 'normal', 'phase_only', or 'absorption_only'. If 'absorption_only', the delta part of the phase of the object will be forced to be 0 after each update. Vice versa for 'phase_only'. |
non_negativity |
Bool | False |
Whether to enforce non-negative constraint. Useful only when unknown_type is delta_beta. |
Forward model¶
| Arg name | Type | Default | Description |
|---|---|---|---|
forward_model |
'auto' or adorym.ForwardModel class |
'auto' |
Forward model class. Use 'auto' to let the program automatically determine forward model from other parameters. |
forward_algorithm |
String | 'fresnel'' |
Choose from 'fresnel' and 'ctf'. |
ctf_lg_kappa |
Float | 1.7 | The natural log of the proportional coefficient between delta and beta, i.e., kappa = 10 ** ctf_lg_kappa; beta_slice = delta_slice * kappa. Only useful when optimize_ctf_lg_kappa is True, in which case the object will be constrained to be homogeneous. Otherwise, delta and beta are reconstructed independently and this argument is ignored. |
binning |
Int | 1 | The number of axial slices to be binned (i.e., to be treated as line integrals) during multislice propagation. |
pure_projection |
Bool | False |
Set to True to model the propagation through the entire object as a simple line projection, not using multislice at all. |
two_d_mode |
Bool | False |
If the HDF5 dataset contains multiple viewing angles (i.e., the length of the first dimension is larger than 1), setting two_d_mode to True will let the program to treat it as a single-angle dataset, with the only angle being the first one. Set to True automatically if the last dimension of obj_size is 1. |
probe_type |
String | 'gaussian' |
Choose from 'gaussian', 'plane', 'ifft', 'aperture_defocus', and 'supplied'. The method of initializing the probe function. Some options requires additional inputs from user. For more details, see table below. |
probe_extra_defocus_cm |
Float | None |
If not None, the probe will be defocused further by the specified distance in cm. |
n_probe_modes |
Int | 1 | Number of probe modes. |
rescale_probe_intensity |
Bool | True |
Scale the probe function so that its integrated power spectrum (related to the total number of photons) matches that of the raw data. |
loss_function_type |
String | 'lsq' |
Choose from 'lsq' or 'poisson'. Whether to use a least square term or a Poisson maximum likelihood term to measure the mismatch of predicted intensity. |
poisson_multiplier |
Float | 1 | Intensity scaling factor in Poisson loss function. If intensity data is normalized, this should be the average number of incident photons per pixel. |
safe_zone_width |
Int | None |
If not None, the object and probe tiles will be enlarged (through either selecting a larger area or padding) before propagation, and the enlarged parts are discarded after propagation. |
scale_ri_by_k |
Bool | True |
Whether to add in the factor k = 2*pi/lambda when evaluating exp(-iknz). Setting this argument to False may help fix numnerical instability problems. |
sign_convention |
Int | 1 | Choose from 1 and -1. Determines whether to use the exp(ikz) convention or exp(-ikz) convention. The reconstructed phase in these two cases will be numerically inverted to each other. |
| Value of ``probe_type`` | Options | Description |
|---|---|---|
'gaussian' |
probe_mag_sigma, probe_phase_sigma, probe_phase_max |
Initialize with a Gaussian probe. The Gaussian spreads, or the *sigma values, are in pixel. Magnitude max is 1 by default. |
'aperture_defocus' |
aperture_radius, beamstop_radius, probe_defocus_cm |
Initialize the probe by defocuing an aperture function. All radii are in pixels (on the object frame). A circular aperture (if beamstop_radius == 0) or a ring aperture (if 0 < beamstop_radius < aperture_radius) is generated and then Fresnel defocused to created the initial probe. |
'ifft' |
Initialize the probe by taking the average of all diffraction patterns, performing an IFFT, and take the moduli. | |
'supplied' |
probe_initial |
Provide a List of Arrays: [probe_mag, probe_phase]. If there are multiple probe modes, each of the arrays should be of shape [n_probe_modes, len_probe_y, len_probe_x]. |
I/O¶
| Arg name | Type | Default | Description |
|---|---|---|---|
save_path |
String | '.' |
Directory that contains the raw data HDF5. If it is in the same folder as the execution script, put '.'. |
output_folder |
String | None |
Name of the folder to place output data. The folder will be assumed to be under save_path, i.e., the actual output directory will be <save_path>/<output_folder>. If None, the folder name will be automatically generated. |
save_intermediate |
Bool | False |
Whether to save the intermediate object (and probe when optimize_probe is True) after each minibatch. |
save_history |
Bool | False |
Useful only if save_intermediate is on, If True, the intermediate output will be saved with a different file name characterized by the current epoch and minibatch number. Otherwise, the intermediate output will be overwritten. |
store_checkpoint |
Bool | True |
Whether to save a checkpoint of the optimizable variables before each minibatch. |
use_checkpoint |
Bool | True |
If set to True, the program initializes the object and/or probe using the checkpoint stored in previous runs. If False or if checkpoint file is not found, start the reconstruction from scratch. |
force_to_use_checkpoint |
Bool | False |
If set to True, when previous checkpoint does not exist or is incomplete, the program raises an error instead of starting from scratch. |
n_batch_per_checkpoint |
Int | 10 | For every how many minibatches should the checkpoint be updated. Large object functions may cause long writing overhead so a larger setting is preferred. |
save_stdout |
Bool | False |
Set to True to save the output messages as a text file. |
Performance¶
| Arg name | Type | Default | Description |
|---|---|---|---|
cpu_only |
Boolean | False |
Set to False to enable GPU. This option is ineffective when backend is autograd. |
gpu_index |
Int | 0 | Index of GPU to use. To use multiple GPUs with multiple MPI ranks, make sure each rank is assigned with a different GPU. |
n_dp_batch |
Int | 20 | Number of tiles to be propagated each time. Values larger than minibatch_size make no difference from setting it equal to minibatch_size. |
distribution_mode |
String or None |
None | Choose from None, 'distributed_object', and 'shared_file', which respectively correspond to data parallel mode, distributed object mode, and H5-mediated low-memory mode. Using the low-memory node requires H5Py built against MPIO-enabled HDF5. |
dist_mode_n_batch_per_update |
Int or None |
None | Update frequency when using distributed object mode. If None, object is updated only after all DPs on an angle are processed. |
precalculate_rotation_coords |
Bool | True |
Whether to calculate rotation transformation coordinates and save them on the hard drive, or calculate them on-the-fly. |
rotate_out_of_loop |
Bool | False |
Applies to simple data parallelism mode only. If True, DP will do rotation outside the loss function and the rotated object function is sent for differentiation. May reduce the number of rotation operations if minibatch_size < n_tiles_per_angle, but object can be updated once only after all tiles on an angle are processed. Also this will save the object-sized gradient array in GPU memory or RAM depending on current device setting. |
Other (non-object) optimizers¶
| Arg name | Type | Default | Description |
|---|---|---|---|
optimize_probe |
Bool | False |
Whether to optimize the probe function. |
probe_learning_rate |
Float | 1e-5 |
Probe optimization step size. |
optimizer_probe |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
optimize_probe_defocusing |
Bool | False |
Whether to optimize the defocusing distance of the probe. |
probe_defocusing_learning_rate |
Float | 1e-5 |
Probe defocusing optimization step size. |
optimizer_probe_defocusing |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
optimize_probe_pos_offset |
Bool | False |
Whether to optimize the offset to probe positions. This is intended to correct for the x-y drifting of the sample stage at different angles. When turned on, the program creates an array with shape [n_rotation_angles, 2]. When processing data from a certain viewing angle, the positions of all diffraction spots are shifted by the value corresponding to that angle. The offset array is optimized by the optimizer along with the object function. |
probe_pos_offset_learning_rate |
Float | 1e-2 |
Probe offset optimization step size. |
optimizer_probe_pos_offset |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
optimize_prj_pos_offset |
Bool | False |
Whether to optimize the offset to projection positions. Note that this is different from optimize_probe_pos_offset: shifting of projections happens after the probe has been modulated by the object function. Useful for refining plane-wave tomography alignment. |
prj_pos_offset_learning_rate |
Float | 1e-2 |
Projection offset optimization step size. |
optimizer_prj_pos_offset |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
optimize_all_probe_pos |
Bool | False |
Whether to optimize the probe positions at all angles. When turned on, the optimizer tries to optimize an array with shape [n_rotation_angles, n_diffraction_spots, 2], which stores the correction values applied to each probe position at all viewing angles. Not recommended for ptychotomography with many viewing angles as it significantly increases the unknwon space to be searched, making the problem less well constrained. |
all_probe_pos_learning_rate |
Float | 1e-2 |
All probe position optimization step size. |
optimizer_all_probe_pos |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
optimize_slice_pos |
Bool | False |
Whether to optimize slice positions. Used for sparse multislice ptychography where slice spacings are not uniform. |
slice_pos_learning_rate |
Float | 1e-4 |
Slice position optimization step size. |
optimizer_slice_pos |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
optimize_free_prop |
Bool | False |
Whether to optimize free propagation distances. |
free_prop_learning_rate |
Float | 1e-2 |
Free propagation distance optimization step size. |
optimizer_free_prop |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
optimize_prj_affine |
Bool | False |
Whether to optimize the affine alignment of holograms. Used for multi-distance holography. |
prj_affine_learning_rate |
Float | 1e-3 |
Affine alignment step size. |
optimizer_prj_affine |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
optimize_tilt |
Bool | False |
Whether to optimize object tilt in all 3 axes. Works only with data parallelism mode. |
tilt_learning_rate |
Float | 1e-3 |
Tilt optimization step size. |
optimizer_tilt |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
initial_tilt |
ndarray |
None |
Initial 3D tilts with shape [3, n_theta]. If not None, 3D tilt will be applied to the object when DP mode is used, regardless whether optimizer_tilt is on or not. |
optimize_ctf_lg_kappa |
Bool | False |
Whether to enable homogeneity constraint and optimize coefficient kappa, where beta_slice = delta_slice * kappa. |
ctf_lg_kappa_learning_rate |
Float | 1e-3 |
kappa optimization step size. |
optimizer_ctf_lg_kappa |
adorym.Optimizer |
None |
Pre-declared optimizer class. If None, a default optimizer will be declared using provided step size and other default parameters. |
other_params_update_delay |
Int | 0 | If larger than 0, updates of above parameters will not happen until the specified number of minibatches are finished. This setting does not apply to object function. |
Other settings¶
| Arg name | Type | Default | Description |
|---|---|---|---|
dynamic_rate |
Bool | True |
Whether to adaptively reduce step size when using GD optimizer. |
debug |
Bool | False |
Whether to enable debugging messages. |
t_max_min |
Float or None |
None | At the end of a batch, terminate the program with s tatus 0 if total time exceeds the set value. Useful for working with supercomputers’ job dependency system, where the dependent may start only if the parent job exits with status 0. |
Forward models¶
Forward models are defined as ForwardModel classes in forward_model.py. All forward model
classes inherent from the parent ForwardModel class:
-
class
adorym.forward_model.ForwardModel(loss_function_type='lsq', distribution_mode=None, device=None, common_vars_dict=None, raw_data_type='magnitude', simulation_mode=False)¶ Bases:
objectThe parent forward model class.
Parameters: - loss_function_type – String. Can be
'lsq'or'poisson'. Type of loss function. - distribution_mode –
None,'distributed_object', or'shared_file'. Mode of parallelization. - device –
Noneor device object. UseNonefor computations on CPU only. - common_vars_dict – Dict. A dictionary of variables that are static (i.e., won’t be optimized and won’t
change with different minibatches). When the ForwardModel class is declared in
reconstruct_ptycho, it can simply be the local namespace of the main module, i.e.,locals(). - raw_data_type – String. Can be
'magnitude'or'intensity'. Type of raw data in HDF5. For line-projection tomography reconstruction where raw data have already been minus-logged, always use'magnitude'.
-
add_regularizers(reg_list)¶
-
get_argument_index(arg)¶
-
get_data(this_i_theta, this_ind_batch, theta_downsample=None, ds_level=1)¶
-
get_loss_function()¶ Overwrite this function for each sub-class, but make sure it contains the explicitly listed arguments. Do not keep differentiable arguments as **kwargs since it will not be recognized by Autograd.
-
get_mismatch_loss(this_pred_batch, this_prj_batch)¶
-
get_regularization_value(obj, device=None)¶
-
loss(this_pred_batch, this_prj_batch, obj)¶ The last-layer loss function that computes (regularized) loss from predicted magnitude. This function is combined with predict to give the full loss function, as done in get_loss_function. Argument this_pred_batch is assumed to be detected magnitude (square root of intensity). If this is not the case for your specific ForwardModel class, override this method.
-
predict(*args, **kwargs)¶ Override and expand this method for each specific class. Generally the method should return the detected magnitude. If this is not the case, you should also override method loss (the last-layer loss) or get_loss_function (the full loss function constructor).
-
update_loss_args(kwargs)¶
- loss_function_type – String. Can be
Currently, the following forward models are built in with Adorym:
-
class
adorym.forward_model.PtychographyModel(loss_function_type='lsq', distribution_mode=None, device=None, common_vars_dict=None, raw_data_type='magnitude', simulation_mode=False)¶ Bases:
adorym.forward_model.ForwardModel-
get_loss_function()¶ Overwrite this function for each sub-class, but make sure it contains the explicitly listed arguments. Do not keep differentiable arguments as **kwargs since it will not be recognized by Autograd.
-
predict(obj, probe_real, probe_imag, probe_defocus_mm, probe_pos_offset, this_i_theta, this_pos_batch, prj, probe_pos_correction, this_ind_batch, tilt_ls, prj_pos_offset)¶ Calculated predicted measurment. The signature of this function exactly matches the function returned by the
get_loss_functionmethod.Parameters: - obj – Array with shape [obj_size_y, obj_size_x, obj_size_z, 2]. Object function.
- probe_real – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Preal part of the probe.
- probe_imag – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Imaginary part of the probe.
- probe_defocus_mm – Probe defocus in mm.
- probe_pos_offset – Array with shape [n_theta, 2]. Probe position offset for each angle, in pixels.
- this_i_theta – Int. Current angle index.
- this_pos_batch – Array with shape [minibatch_size, 2]. Current batch of probe positions.
- prj – HDF5 measurement dataset handler.
- probe_pos_correction – Array with shape [n_probes, 2]. Additive correction values of probe positions.
- this_ind_batch – List of Int. Current batch of tile indices.
- tilt_ls – Array with shape [3, n_theta]. 3-axis tilt to be applied before propagation.
- prj_pos_offset – Array with shape [n_theta, 2]. Projection position offset for each angle, in pixels. Note that this is different from probe_pos_offset - shifting is applied after passing through the object.
Returns: Array with shape [minibatch_size, len_probe_y, len_probe_x]. Magnitude of detected wavefields.
-
-
class
adorym.forward_model.SingleBatchFullfieldModel(loss_function_type='lsq', distribution_mode=None, device=None, common_vars_dict=None, raw_data_type='magnitude', simulation_mode=False)¶ Bases:
adorym.forward_model.PtychographyModel-
predict(obj, probe_real, probe_imag, probe_defocus_mm, probe_pos_offset, this_i_theta, this_pos_batch, prj, probe_pos_correction, this_ind_batch, tilt_ls, prj_pos_offset)¶ Calculated predicted measurment. The signature of this function exactly matches the function returned by the
get_loss_functionmethod.Parameters: - obj – Array with shape [obj_size_y, obj_size_x, obj_size_z, 2]. Object function.
- probe_real – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Preal part of the probe.
- probe_imag – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Imaginary part of the probe.
- probe_defocus_mm – Probe defocus in mm.
- probe_pos_offset – Array with shape [n_theta, 2]. Probe position offset for each angle, in pixels.
- this_i_theta – Int. Current angle index.
- this_pos_batch – Array with shape [minibatch_size, 2]. Current batch of probe positions.
- prj – HDF5 measurement dataset handler.
- probe_pos_correction – Array with shape [n_probes, 2]. Additive correction values of probe positions.
- this_ind_batch – List of Int. Current batch of tile indices.
- tilt_ls – Array with shape [3, n_theta]. 3-axis tilt to be applied before propagation.
Returns: Array with shape [minibatch_size, len_probe_y, len_probe_x]. Magnitude of detected wavefields.
-
-
class
adorym.forward_model.SingleBatchPtychographyModel(loss_function_type='lsq', distribution_mode=None, device=None, common_vars_dict=None, raw_data_type='magnitude', simulation_mode=False)¶ Bases:
adorym.forward_model.PtychographyModel-
predict(obj, probe_real, probe_imag, probe_defocus_mm, probe_pos_offset, this_i_theta, this_pos_batch, prj, probe_pos_correction, this_ind_batch, tilt_ls, prj_pos_offset)¶ Calculated predicted measurment. The signature of this function exactly matches the function returned by the
get_loss_functionmethod.Parameters: - obj – Array with shape [obj_size_y, obj_size_x, obj_size_z, 2]. Object function.
- probe_real – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Preal part of the probe.
- probe_imag – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Imaginary part of the probe.
- probe_defocus_mm – Probe defocus in mm.
- probe_pos_offset – Array with shape [n_theta, 2]. Probe position offset for each angle, in pixels.
- this_i_theta – Int. Current angle index.
- this_pos_batch – Array with shape [minibatch_size, 2]. Current batch of probe positions.
- prj – HDF5 measurement dataset handler.
- probe_pos_correction – Array with shape [n_probes, 2]. Additive correction values of probe positions.
- this_ind_batch – List of Int. Current batch of tile indices.
- tilt_ls – Array with shape [3, n_theta]. 3-axis tilt to be applied before propagation.
Returns: Array with shape [minibatch_size, len_probe_y, len_probe_x]. Magnitude of detected wavefields.
-
-
class
adorym.forward_model.SparseMultisliceModel(loss_function_type='lsq', distribution_mode=None, device=None, common_vars_dict=None, raw_data_type='magnitude', simulation_mode=False)¶ Bases:
adorym.forward_model.ForwardModel-
get_loss_function()¶ Overwrite this function for each sub-class, but make sure it contains the explicitly listed arguments. Do not keep differentiable arguments as **kwargs since it will not be recognized by Autograd.
-
predict(obj, probe_real, probe_imag, probe_defocus_mm, probe_pos_offset, this_i_theta, this_pos_batch, prj, probe_pos_correction, this_ind_batch, slice_pos_cm_ls, prj_pos_offset)¶ Calculated predicted measurment. The signature of this function exactly matches the function returned by the
get_loss_functionmethod.Parameters: - obj – Array with shape [obj_size_y, obj_size_x, obj_size_z, 2]. Object function.
- probe_real – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Preal part of the probe.
- probe_imag – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Imaginary part of the probe.
- probe_defocus_mm – Probe defocus in mm.
- probe_pos_offset – Array with shape [n_theta, 2]. Probe position offset for each angle, in pixels.
- this_i_theta – Int. Current angle index.
- this_pos_batch – Array with shape [minibatch_size, 2]. Current batch of probe positions.
- prj – HDF5 measurement dataset handler.
- probe_pos_correction – Array with shape [n_probes, 2]. Additive correction values of probe positions.
- this_ind_batch – List of Int. Current batch of tile indices.
- slice_pos_cm_ls – Array of Float. Positions of object slices in cm.
- prj_pos_offset – Array with shape [n_theta, 2]. Projection position offset for each angle, in pixels. Note that this is different from probe_pos_offset - shifting is applied after passing through the object.
Returns: Array with shape [minibatch_size, len_probe_y, len_probe_x]. Magnitude of detected wavefields.
-
-
class
adorym.forward_model.MultiDistModel(loss_function_type='lsq', distribution_mode=None, device=None, common_vars_dict=None, raw_data_type='magnitude', simulation_mode=False)¶ Bases:
adorym.forward_model.ForwardModel-
get_loss_function()¶ Overwrite this function for each sub-class, but make sure it contains the explicitly listed arguments. Do not keep differentiable arguments as **kwargs since it will not be recognized by Autograd.
-
predict(obj, probe_real, probe_imag, probe_defocus_mm, probe_pos_offset, this_i_theta, this_pos_batch, prj, probe_pos_correction, this_ind_batch, free_prop_cm, safe_zone_width, prj_affine_ls, ctf_lg_kappa, prj_pos_offset)¶ Calculated predicted measurment. The signature of this function exactly matches the function returned by the
get_loss_functionmethod.Parameters: - obj – Array with shape [obj_size_y, obj_size_x, obj_size_z, 2]. Object function.
- probe_real – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Preal part of the probe.
- probe_imag – Array with shape [n_probe_modes, len_probe_y, len_probe_x]. Imaginary part of the probe.
- probe_defocus_mm – Probe defocus in mm.
- probe_pos_offset – Array with shape [n_theta, 2]. Probe position offset for each angle, in pixels.
- this_i_theta – Int. Current angle index.
- this_pos_batch – Array with shape [minibatch_size, 2]. Current batch of probe positions.
- prj – HDF5 measurement dataset handler.
- probe_pos_correction – Array with shape [n_probes, 2]. Additive correction values of probe positions.
- this_ind_batch – List of Int. Current batch of tile indices.
- free_prop_cm – Array of Float. Propagation distances for all holograms.
- safe_zone_width – Int. Width of safe zone to be padded to probes before propagation to prevent fringe wrapping.
- prj_affine_ls – Array with shape [n_dists, 2, 3]. Affine transform matrices for each hologram.
- ctf_lg_kappa – Log10 of kappa, the relation coefficient between delta and beta. If homogeneous assumption is not used, set as None.
- prj_pos_offset – Array with shape [n_theta, 2]. Projection position offset for each angle, in pixels. Note that this is different from probe_pos_offset - shifting is applied after passing through the object.
Returns: Array with shape [minibatch_size, len_probe_y, len_probe_x]. Magnitude of detected wavefields.
-
For all forward models, this_i_theta, this_pos_batch, and this_ind_batch
are essential since they inform the program the right data to read for the current minibatch.
Some static variables, like the basic probe positions (but not the correction values to
probe positions), are not included as an argument, but are directly fetched from common_vars_dict.
Regularizers¶
Regularizer classes can be declared and passed as a list to reconstruct_ptychography via
the regularizers argument. The list of regularizers is binded to the ForwardModel class
through forward_model.add_regularizers.
-
class
adorym.regularizers.Regularizer(unknown_type='delta_beta')¶ Bases:
objectParent regularizer class.
Parameters: - unknown_type – String. Can be
'delta_beta'or'real_imag'. - device – Device object or
None.
-
get_value(obj, device=None, **kwargs)¶
- unknown_type – String. Can be
The following regularizers are available:
-
class
adorym.regularizers.L1Regularizer(alpha_d, alpha_b, unknown_type='delta_beta')¶ Bases:
adorym.regularizers.RegularizerL1-norm regularizer.
Parameters: - alpha_d – Weight of l1-norm of delta or real part.
- alpha_b – Weight of l1-norm of beta or imaginary part.
-
get_value(obj, device=None, **kwargs)¶
-
class
adorym.regularizers.ReweightedL1Regularizer(alpha_d, alpha_b, unknown_type='delta_beta')¶ Bases:
adorym.regularizers.RegularizerReweighted l1-norm regularizer.
Parameters: - alpha_d – Weight of l1-norm of delta or real part.
- alpha_b – Weight of l1-norm of beta or imaginary part.
-
get_value(obj, device=None, **kwargs)¶
-
update_l1_weight(weight_l1)¶
-
class
adorym.regularizers.TVRegularizer(gamma, unknown_type='delta_beta')¶ Bases:
adorym.regularizers.RegularizerTotal variation regularizer.
Parameters: gamma – Weight of TV term. -
get_value(obj, distribution_mode=None, device=None, **kwargs)¶
-
An alternative to explicitly declaring the regularizer objects is to supply the alpha_d,
alpha_b, gamma, and reweighted_l1 arguments in reconstruct_ptychography.
Optimizers¶
When setting the optimizer for the object function, users can provide
the name of the optimizer (see
Object optimizer options) and the step
size of that parameter as the only hyperparameter.
For other refinable parameters, users may use the default optimizer
type, only specifying the step size. This can
be limited when one wants to try different types of optimizers for
non-object variables or to tune optimizer hyperparameters
other than the step size. Therefore, you may also explicitly declare
the optimizer, and pass the adorym.Optimizer
object to the main fucntion reconstruct_ptychography.
For now, ScipyOptimizer can only be used for the object function.
Below is the API reference of the general Optimizer class:
-
class
adorym.optimizers.Optimizer(name, output_folder='.', params_list=(), distribution_mode=None, options_dict=None, forward_model=None)¶ Bases:
objectParent class of optimizers.
Parameters: - name – Name of the optimizer. It is currently used to (1) match the optimizer to special handling rules defined in optimizers.update_parameters, optimizers.update_parameter_gradients, optimizers.create_parameter_output_folders, and optimizers.output_intermediate_parameters, and (2) to locate the optimized variable in predict function’s argument list in ScipyOptmizer. If the optimizer is created for preset variables (e.g., probe_pos_correction), the name can be any arbitrary string since Adorym will forcefully set the names to the default names for these variables. If the optimizer if created for user-defined optimizable parameters, make sure the name is the same as the name of the variable being optimized, and matches the rules defined in the aforementioned functions, if any.
- output_folder – String. Path to the output folder. This should be the combination of save_path and output_folder passed to reconstruct_ptytchography. This path will be the location to save/read checkpoints of optimizer parameters.
- distribution_mode – None or String. Should match the value passed to reconstruct_ptychography.
- options_dict – Dict. A dictionary of optimizer hyperparameters. The options differ depending on the type of optimizers. See table below for a thorough reference.
- params_list – List of str; a list of optimizer parameters provided in strings. (Only used for parent class. Child Optimizer classes do NOT take this argument.)
-
convert_gradient(gradient)¶
-
create_container(whole_object_size, use_checkpoint, device_obj, use_numpy=False, dtype='float32')¶ Parameters: whole_object_size – List of int; 4-D vector for object function (including 2 channels), or a 3-D vector for probe, or a 1-D scalar for other variables. Channel must be the last domension. Parameter arrays will be created following exactly whole_object_size.
-
create_distributed_param_arrays(whole_object_size, use_numpy=False, dtype='float32')¶
-
create_file_objects(whole_object_size, use_checkpoint=False)¶
-
create_param_arrays(whole_object_size, device=None, use_numpy=False)¶
-
get_array_slicer(slicer)¶
-
get_params_from_file(this_pos_batch=None, probe_size=None)¶
-
read_chunks_from_distributed_object(probe_pos, this_ind_batch_allranks, minibatch_size, probe_size, device=None, unknown_type='delta_beta', apply_to_arr_rot=False, dtype='float32', n_split='auto')¶
-
restore_distributed_param_arrays_from_checkpoint(device=None, use_numpy=False, dtype='float32')¶
-
restore_param_arrays_from_checkpoint(device=None, use_numpy=False)¶
-
rotate_arrays(coords, interpolation='bilinear', overwrite_arr=False)¶
-
rotate_files(coords, interpolation='bilinear')¶
-
save_distributed_param_arrays_to_checkpoint(use_numpy=True)¶
-
save_param_arrays_to_checkpoint(use_numpy=False)¶
-
set_index_in_grad_return(ind)¶
-
sync_chunks_to_distributed_object(arr, probe_pos, this_ind_batch_allranks, minibatch_size, probe_size, dtype='float32', n_split='auto')¶
-
write_params_to_file(this_pos_batch=None, probe_size=None, n_ranks=1)¶
The acceptable values of options_dict varies with different optimizer types.
The table below lists the available options:
| Optimizer | ``options_dict`` and default values |
|---|---|
GDOptimizer |
step_size=0.001, dynamic_rate=True, first_downrate_iteration=92 |
AdamOptimizer |
step_size=0.001, b1=0.9, b2=0.999, eps=1e-7 |
MomentumOptimizer |
step_size=0.001, gamma=0.9 |
CurveballOptimizer |
alpha=1.0 |
CGOptimizer |
step_size=1.0, linesearch_type='adaptive', max_backtracking_iter=None |
ScipyOptimizer* |
step_size=1.e2, method='CG', options=None** |
* ScipyOptimizer needs Hessian-vector product when method is one
of Newton-CG, trust-ncg, trust-krylov,
and trust-constr. In these cases, the HVP is approximated using
Gauss-Newton method.
** For valid values of method and options, refer to the
documentation of scipy.optimize.minimize.
Generally, the apply_gradient method in each optimizer class is called to update
the variable. Some Special optimizers, like CurveballOptimizer, may require additional
calls to other methods in order to compute intermediate parameters.
-
class
adorym.optimizers.CurveballOptimizer(name, output_folder='.', distribution_mode=None, options_dict=None, forward_model=None)¶ Bases:
adorym.optimizers.OptimizerGauss-Newton second-order optimizer implemented using the Curveball algorithm: Henriques, J. F., Ehrhardt, S., Albanie, S. & Vedaldi, A. Small steps and giant leaps: Minimal Newton solvers for Deep Learning. arXiv (2018). This code is adapted from https://github.com/saugatkandel/sopt. When working with DO, dz is synchronized in place of gradient.
-
apply_gradient(x, gradient, i_batch, alpha=1, use_numpy=False, params_slicer=None)¶ Use calculated gradient to update the variable being optimized. :param x: Array or Tensor of the optimized variable. :param gradient: Array or adorym.Gradient. If optimizer is CG, the ForwardModel instance (which is needed for
providing loss function for line search) can be supplied through the Gradient instance. Otherwise, it must be specified when the optimizer is instantiated.Parameters: i_batch – Int. User-specifiable step number. When minibatching localized data using optimizers like Adam, i_batch may be preferably up-counted only when all voxels of the object are updated with non-zero gradient.
-
calculate_beta_rho(differentiator, use_numpy=False)¶ Parameters are calculated using chunks when working with DO. In DP mode, self.dz_chunk and self.z_chunk should match object size.
-
calculate_dz(differentiator, use_numpy=False)¶ In DO, dz will be synchronized as Gradient class after this step.
-
calculate_update_vector(dz, ss)¶ In DO, this is done for slabs.
-
update_lambda(forward_model, forward_args)¶
-
Output¶
During runtime, Adorym may create a folder named
arrsize_?_?_?_ntheta_? in the current working directory, which
saves
the precalculated coordinates for rotation transformation. Other than
that, all outputs will be written in <save_path>/<output_folder>,
which is organized as shown in the chart below:
output_folder
|___ convergence
| |___ loss_rank_0.txt // Record of the loss value after
| |___ loss_rank_1.txt // each update coming from each process.
| |___ ...
|___ intermediate
| |___ object
| | |___ obj_mag(delta)_0_0.tiff
| | |___ obj_phase(beta)_0_0.tiff
| | |___ ...
| |___ probe
| | |___ probe_mag_0_0.tiff
| | |___ probe_phase_0_0.tiff
| | |___ ...
| |___ probe_pos (if optimize_all_probe_pos is True)
| | |___ probe_pos_correction_0_0_0.txt
| | |___ ...
| ...
|___ obj_delta_ds_1.tiff (or obj_mag_ds_1.tiff)
|___ obj_beta_ds_1.tiff (or obj_phase_ds_1.tiff)
|___ probe_mag_ds_1.tiff
|___ probe_phase_ds_1.tiff
|___ summary.txt // Summary of parameter settings.
|___ checkpoint.txt // Exists if store_checkpoint is True.
|___ obj_checkpoint.npy // Exists if store_checkpoint is True.
|___ opt_params_checkpoint.npy // Exists if store_checkpoint is True and optimizer has parameters.
By default, all image outputs are in 32-bit floating points which can be opened and viewed with ImageJ.