aspire.abinitio package¶
Submodules¶
aspire.abinitio.commonline_base module¶
- class aspire.abinitio.commonline_base.CLOrient3D(src, n_rad=None, n_theta=360, n_check=None, hist_bin_width=3, full_width=6, max_shift=0.15, shift_step=1, offsets_max_shift=None, offsets_shift_step=None, mask=True)¶
- Bases: - object- Define a base class for estimating 3D orientations using common lines methods - Initialize an object for estimating 3D orientations using common lines. - Parameters:
- src – The source object of 2D denoised or class-averaged images. 
- n_rad – The number of points in the radial direction. If None, n_rad will default to the ceiling of half the resolution of the source. 
- n_theta – The number of points in the theta direction. This value must be even. Default is 360. 
- n_check – For each image/projection find its common-lines with n_check images. If n_check is less than the total number of images, a random subset of n_check images is used. 
- max_shift – Determines maximum range for shifts for common-line detection as a proportion of the resolution. Default is 0.15. 
- shift_step – Resolution of shift estimation for common-line detection in pixels. Default is 1 pixel. 
- offsets_max_shift – Determines maximum range for shifts for 2D offset estimation as a proportion of the resolution. Default None inherits from max_shift. 
- offsets_shift_step – Resolution of shift estimation for 2D offset estimation in pixels. Default None inherits from shift_step. 
- hist_bin_width – Bin width in smoothing histogram (degrees). 
- full_width – Selection width around smoothed histogram peak (degrees). adaptive will attempt to automatically find the smallest number of `hist_bin_width`s required to find at least one valid image index. 
- mask – Option to mask src.images with a fuzzy mask (boolean). Default, True, applies a mask. 
 
 - build_clmatrix()¶
- Build common-lines matrix from Fourier stack of 2D images - Wrapper for cpu/gpu dispatch. 
 - build_clmatrix_cu()¶
- Build common-lines matrix from Fourier stack of 2D images 
 - build_clmatrix_host()¶
- Build common-lines matrix from Fourier stack of 2D images 
 - property clmatrix¶
- Returns Common Lines Matrix. - Computes if clmatrix is None. - Returns:
- Common Lines Matrix 
 
 - estimate(**kwargs)¶
- Estimate orientation and shifts for all 2D images. - Returns:
- (rotations, shifts) 
 
 - estimate_rotations()¶
- Estimate orientation matrices for all 2D images - Subclasses should implement this function. 
 - estimate_shifts(equations_factor=1, max_memory=4000)¶
- Estimate 2D shifts in images - This function computes 2D shifts in x, y of images by solving the least-squares equations to Ax = b. A on the left-hand side is a sparse matrix representing precomputed coefficients of shift equations; and on the right-side, b is estimated 1D shifts along the theta direction between two Fourier rays (one in image i and the other in image j). Each row of shift equations contains four unknowns, shifts in x, y for a pair of images. The detailed implementation can be found in the book chapter as below: Y. Shkolnisky and A. Singer, Center of Mass Operators for Cryo-EM - Theory and Implementation, Modeling Nanoscale Imaging in Electron Microscopy, T. Vogt, W. Dahmen, and P. Binev (Eds.) Nanostructure Science and Technology Series, Springer, 2012, pp. 147–177 - Parameters:
- equations_factor – The factor to rescale the number of shift equations (=1 in default) 
- max_memory – If there are N images and N_check selected to check for common lines, then the exact system of equations solved for the shifts is of size 2N x N(N_check-1)/2 (2N unknowns and N(N_check-1)/2 equations). This may be too big if N is large. The algorithm will use equations_factor times the total number of equations if the resulting total number of memory requirements is less than max_memory (in megabytes); otherwise it will reduce the number of equations by approximation to fit in max_memory. 
 
 
 - property pf¶
 - property rotations¶
- Returns estimated rotations. - Computes if rotations is None. - Returns:
- Estimated rotations 
 
 - property shifts¶
- Returns estimated shifts. - Computes if shifts is None. - Returns:
- Estimated shifts 
 
 
aspire.abinitio.commonline_c2 module¶
- class aspire.abinitio.commonline_c2.CLSymmetryC2(src, n_rad=None, n_theta=360, max_shift=0.15, shift_step=1, epsilon=0.001, max_iters=1000, degree_res=1, min_dist_cls=25, seed=None, mask=True)¶
- Bases: - CLSymmetryC3C4- Define a class to estimate 3D orientations using common lines methods for molecules with C2 cyclic symmetry. - The related publications are listed below: X. Cheng, Random Matrices in High-dimensional Data Analysis, PhD thesis, Princeton University, (2013). - G. Pragier and Y. Shkolnisky, A Common Lines Approach for Abinitio Modeling of Cyclically Symmetric Molecules, Inverse Problems, 35, 124005, (2019). - Y. Shkolnisky, and A. Singer, Viewing Direction Estimation in Cryo-EM Using Synchronization, SIAM J. Imaging Sciences, 5, 1088-1110 (2012). - A. Singer, R. R. Coifman, F. J. Sigworth, D. W. Chester, Y. Shkolnisky, Detecting Consistent Common Lines in Cryo-EM by Voting, Journal of Structural Biology, 169, 312-322 (2010). - Initialize object for estimating 3D orientations for molecules with C2 symmetry. - Parameters:
- src – The source object of 2D denoised or class-averaged images with metadata 
- n_rad – The number of points in the radial direction 
- n_theta – The number of points in the theta direction. Default = 360. 
- max_shift – Maximum range for shifts as a proportion of resolution. Default = 0.15. 
- shift_step – Resolution of shift estimation in pixels. Default = 1 pixel. 
- epsilon – Tolerance for the power method. 
- max_iter – Maximum iterations for the power method. 
- degree_res – Degree resolution for estimating in-plane rotations. 
- min_dist_cls – Minimum distance between mutual common-lines. Default = 25 degrees. 
- seed – Optional seed for RNG. 
- mask – Option to mask src.images with a fuzzy mask (boolean). Default, True, applies a mask. 
 
 - build_clmatrix()¶
- Build common-lines matrix for molecules with C2 symmetry from Fourier stack of 2D images. This consists of finding for each pair of images the two common-lines induced by the 2-fold symmetry. 
 - estimate_rotations()¶
- Estimate rotation matrices for molecules with C2 symmetry. 
 
aspire.abinitio.commonline_c3_c4 module¶
- class aspire.abinitio.commonline_c3_c4.CLSymmetryC3C4(src, symmetry=None, n_rad=None, n_theta=360, max_shift=0.15, shift_step=1, epsilon=0.01, max_iters=1000, degree_res=1, seed=None, mask=True)¶
- Bases: - CLOrient3D,- SyncVotingMixin- Define a class to estimate 3D orientations using common lines methods for molecules with C3 and C4 cyclic symmetry. - The related publications are listed below: G. Pragier and Y. Shkolnisky, A Common Lines Approach for Abinitio Modeling of Cyclically Symmetric Molecules, Inverse Problems, 35, 124005, (2019). - Y. Shkolnisky, and A. Singer, Viewing Direction Estimation in Cryo-EM Using Synchronization, SIAM J. Imaging Sciences, 5, 1088-1110 (2012). - A. Singer, R. R. Coifman, F. J. Sigworth, D. W. Chester, Y. Shkolnisky, Detecting Consistent Common Lines in Cryo-EM by Voting, Journal of Structural Biology, 169, 312-322 (2010). - Initialize object for estimating 3D orientations for molecules with C3 and C4 symmetry. - Parameters:
- src – The source object of 2D denoised or class-averaged images with metadata 
- symmetry – A string, ‘C3’ or ‘C4’, indicating the symmetry type. 
- n_rad – The number of points in the radial direction 
- n_theta – The number of points in the theta direction. Default = 360. 
- max_shift – Maximum range for shifts as a proportion of resolution. Default = 0.15. 
- shift_step – Resolution of shift estimation in pixels. Default = 1 pixel. 
- epsilon – Tolerance for the power method. 
- max_iter – Maximum iterations for the power method. 
- degree_res – Degree resolution for estimating in-plane rotations. 
- seed – Optional seed for RNG. 
- mask – Option to mask src.images with a fuzzy mask (boolean). Default, True, applies a mask. 
 
 - static cl_angles_to_ind(cl_angles, n_theta)¶
 - estimate_rotations()¶
- Estimate rotation matrices for molecules with C3 or C4 symmetry. - Returns:
- Array of rotation matrices, size n_imgx3x3. 
 
 - static g_sync(rots, order, rots_gt)¶
- Every estimated rotation might be a version of the ground truth rotation rotated by g^{s_i}, where s_i = 0, 1, …, order. This method synchronizes the ground truth rotations so that only a single global rotation need be applied to all estimates for error analysis. - Parameters:
- rots – Estimated rotation matrices 
- order – The cyclic order asssociated with the symmetry of the underlying molecule. 
- rots_gt – Ground truth rotation matrices. 
 
- Returns:
- g-synchronized ground truth rotations. 
 
 
aspire.abinitio.commonline_cn module¶
- class aspire.abinitio.commonline_cn.CLSymmetryCn(src, symmetry=None, n_rad=None, n_theta=360, max_shift=0.15, shift_step=1, epsilon=0.001, max_iters=1000, degree_res=1, n_points_sphere=500, equator_threshold=10, seed=None, mask=True)¶
- Bases: - CLSymmetryC3C4- Define a class to estimate 3D orientations using common lines methods for molecules with Cn cyclic symmetry, with n>4. - Initialize object for estimating 3D orientations for molecules with Cn symmetry, n>4. - Parameters:
- src – The source object of 2D denoised or class-averaged images with metadata 
- symmetry – A string, ‘Cn’, indicating the symmetry type. 
- n_rad – The number of points in the radial direction. 
- n_theta – The number of points in the theta direction. Default = 360. 
- max_shift – Maximum range for shifts as a proportion of resolution. Default = 0.15. 
- shift_step – Resolution of shift estimation in pixels. Default = 1 pixel. 
- epsilon – Tolerance for the power method. 
- max_iter – Maximum iterations for the power method. 
- degree_res – Degree resolution for estimating in-plane rotations. 
- n_points_sphere – The number of candidate rotations used to estimate viewing directions. 
- equator_threshold – Threshold for removing candidate rotations within equator_threshold degrees of being an equator image. Default is 10 degrees. 
- seed – Optional seed for RNG. 
- mask – Option to mask src.images with a fuzzy mask (boolean). Default, True, applies a mask. 
 
 - estimate_rotations()¶
- Estimate rotation matrices for molecules with Cn symmetry, n > 4. - Returns:
- Array of rotation matrices, size n_imgx3x3. 
 
 - static generate_candidate_rots(n, equator_threshold, order, degree_res, seed)¶
- Generate random rotations that exclude rotations inducing equator images for use as candidates in the CLSymmetryCn algorithm. - Parameters:
- n – Number of rotations to generate. 
- equator_threshold – Angular distance from equator (in degrees). 
- order – Cyclic order of underlying molecule. 
- degree_res – Degree resolution for in-plane rotations. 
- seed – Random seed. 
 
- Returns:
- Candidate rotations, In-plane rotations 
 
 - static relative_rots_to_cl_indices(relative_rots, n_theta)¶
- Given a set of relative rotations between pairs of images produce the common-line indices for each pair. - Parameters:
- relative_rots – The n x 3 x 3 relative rotations between pairs of images. 
- n_theta – The theta resolution for common-line indices. 
 
- Returns:
- Common-line indices c1s, c2s each length n. 
 
 
- class aspire.abinitio.commonline_cn.MeanOuterProductEstimator¶
- Bases: - object- Incrementally accumulate outer product entries of unknown conjugation. - dtype¶
- alias of - float64
 - push(V)¶
- Given V, accumulate entries into a running sum of J-synchronized entries. 
 - synchronized_mean()¶
- Calculate the mean of synchronized outer product estimates. 
 
aspire.abinitio.commonline_d2 module¶
- class aspire.abinitio.commonline_d2.CLSymmetryD2(src, n_rad=None, n_theta=360, max_shift=0.15, shift_step=1, grid_res=1200, inplane_res=5, eq_min_dist=7, epsilon=0.01, seed=None, mask=True)¶
- Bases: - CLOrient3D- Define a class to estimate 3D orientations using common lines methods for molecules with D2 (dihedral) symmetry. - Corresponding publication: E. Rosen and Y. Shkolnisky, Common lines ab-initio reconstruction of D2-symmetric molecules, SIAM Journal on Imaging Sciences, volume 13-4, p. 1898-1994, 2020 - Initialize object for estimating 3D orientations for molecules with D2 symmetry. - Parameters:
- src – The source object of 2D denoised or class-averaged images with metadata 
- n_rad – The number of points in the radial direction of Fourier image. 
- n_theta – The number of points in the theta direction of Fourier image. Default = 360. 
- max_shift – Maximum range for shifts as a proportion of resolution. Default = 0.15. 
- shift_step – Resolution of shift estimation in pixels. Default = 1 pixel. 
- grid_res – Number of sampling points on sphere for projetion directions. These are generated using the Saaf-Kuijlaars algorithm. Default value is 1200. 
- inplane_res – The sampling resolution of in-plane rotations for each projection direction. Default value is 5 degrees. 
- eq_min_dist – Width of strip around equator projection directions from which we do not sample directions. Default value is 7 degrees. 
- epsilon – Tolerance for J-synchronization power method. 
- seed – Optional seed for RNG. 
- mask – Option to mask src.images with a fuzzy mask (boolean). Default, True, applies a mask. 
 
 - estimate_rotations()¶
- Estimate rotation matrices for molecules with D2 symmetry. Sets the attribute self.rotations with an array of estimated rotation matrices, size src.nx3x3. 
 
aspire.abinitio.commonline_ev module¶
- class aspire.abinitio.commonline_ev.CommLineEV(src)¶
- Bases: - CLOrient3D- Class to estimate 3D orientations using Eigenvector method [SS11] - constructor of an object for estimating 3D orientations - estimate()¶
- perform estimation of orientations 
 - output()¶
- Output the 3D orientations 
 
aspire.abinitio.commonline_gcar module¶
- class aspire.abinitio.commonline_gcar.CommLineGCAR(src)¶
- Bases: - CLOrient3D- Define a derived class to estimate 3D orientations using Globally Consistent Angular Reconstitution described as below: R. Coifman, Y. Shkolnisky, F. J. Sigworth, and A. Singer, Reference Free Structure Determination through Eigenvestors of Center of Mass Operators, Applied and Computational Harmonic Analysis, 28, 296-312 (2010). - constructor of an object for estimating 3D oreintations - estimate()¶
- perform estimation of orientations 
 - output()¶
- Output the 3D orientations 
 
aspire.abinitio.commonline_irls module¶
- class aspire.abinitio.commonline_irls.CommonlineIRLS(src, *args, num_itrs=10, eps_weighting=0.001, alpha=None, max_rankZ=None, max_rankW=None, **kwargs)¶
- Bases: - CommonlineLUD- Estimate 3D orientations using Iteratively Reweighted Least Squares (IRLS) as described in the following publication: L. Wang, A. Singer, and Z. Wen, Orientation Determination of Cryo-EM Images Using Least Unsquared Deviations, SIAM J. Imaging Sciences, 6, 2450-2483 (2013). - Initialize a class for estimating 3D orientations using an IRLS-based optimization. See CommonlineLUD for additional arguments. - Parameters:
- num_itrs – Number of iterations for iterative reweighting. Default is 10. 
- eps_weighting – Regularization value for reweighting factor. Default is 1e-3. 
- alpha – Spectral norm constraint for IRLS algorithm. Default is None, which does not apply a spectral norm constraint. To apply a spectral norm constraint provide a value in the range [2/3, 1), 2/3 is recommended. 
- max_rankZ – Maximum rank used for projecting the Z matrix (for adaptive projection). If None, defaults to max(6, n_img // 4). 
- max_rankW – Maximum rank used for projecting the W matrix (for adaptive projection). If None, defaults to max(6, n_img // 4). 
 
 - estimate_rotations()¶
- Estimate rotation matrices using the common lines method with IRLS optimization. 
 
aspire.abinitio.commonline_lud module¶
- class aspire.abinitio.commonline_lud.CommonlineLUD(src, alpha=None, tol=0.001, mu=1, gam=1.618, eps=None, maxit=1000, adp_proj=True, max_rankZ=None, max_rankW=None, adp_mu=True, dec_mu=0.5, inc_mu=2, mu_min=0.0001, mu_max=10000.0, min_mu_itr=5, max_mu_itr=20, delta_mu_l=0.1, delta_mu_u=10, **kwargs)¶
- Bases: - CommonlineSDP- Define a derived class to estimate 3D orientations using Least Unsquared Deviations as described in the following publication: L. Wang, A. Singer, and Z. Wen, Orientation Determination of Cryo-EM Images Using Least Unsquared Deviations, SIAM J. Imaging Sciences, 6, 2450-2483 (2013). - Initialize a class for estimating 3D orientations using a Least Unsquared Deviations algorithm. - This class extends the CLOrient3D class, inheriting its initialization parameters. Additional parameters detailed below. - Parameters:
- alpha – Spectral norm constraint for ADMM algorithm. Default is None, which does not apply a spectral norm constraint. To apply a spectral norm constraint provide a value in the range [2/3, 1), 2/3 is recommended. 
- tol – Tolerance for convergence. The algorithm stops when conditions reach this threshold. Default is 1e-3. 
- mu – The penalty parameter (or dual variable scaling factor) in the optimization problem. Default is 1. 
- gam – Relaxation factor for updating variables in the algorithm (typically between 1 and 2). Default is 1.618. 
- eps – Small positive value used to filter out negligible eigenvalues. Default is 1e-5 if src.dtype is singles, otherwise 1e-12. 
- maxit – Maximum number of iterations allowed for the algorithm. Default is 1000. 
- adp_proj – Flag for using adaptive projection during eigenvalue computation: - True: Adaptive rank selection (Default). - False: Full eigenvalue decomposition. 
- max_rankZ – Maximum rank used for projecting the Z matrix (for adaptive projection). If None, defaults to max(6, n_img // 2). 
- max_rankW – Maximum rank used for projecting the W matrix (for adaptive projection). If None, defaults to max(6, n_img // 2). 
- adp_mu – Adaptive adjustment of the penalty parameter mu: - True: Enabled (Default). - False: Disabled. 
- dec_mu – Scaling factor for decreasing mu. Default is 0.5. 
- inc_mu – Scaling factor for increasing mu. Default is 2. 
- mu_min – Minimum allowable value for mu. Default is 1e-4. 
- mu_max – Maximum allowable value for mu. Default is 1e4. 
- min_mu_itr – Minimum number of iterations before mu is adjusted. Default is 5. 
- max_mu_itr – Maximum number of iterations allowed for mu adjustment. Default is 20. 
- delta_mu_l – Lower bound for relative drop ratio to trigger a decrease in mu. Default is 0.1. 
- delta_mu_u – Upper bound for relative drop ratio to trigger an increase in mu. Default is 10. 
 
 - estimate_rotations()¶
- Estimate rotation matrices using the common lines method with LUD optimization. 
 
aspire.abinitio.commonline_sdp module¶
- class aspire.abinitio.commonline_sdp.CommonlineSDP(src, n_rad=None, n_theta=360, n_check=None, hist_bin_width=3, full_width=6, max_shift=0.15, shift_step=1, offsets_max_shift=None, offsets_shift_step=None, mask=True)¶
- Bases: - CLOrient3D- Class to estimate 3D orientations using semi-definite programming. - See the following publication for more details: - A. Singer and Y. Shkolnisky, “Three-Dimensional Structure Determination from Common Lines in Cryo-EM - by Eigenvectors and Semidefinite Programming” - SIAM J. Imaging Sciences, Vol. 4, No. 2, (2011): 543-572. doi:10.1137/090767777 - Initialize an object for estimating 3D orientations using common lines. - Parameters:
- src – The source object of 2D denoised or class-averaged images. 
- n_rad – The number of points in the radial direction. If None, n_rad will default to the ceiling of half the resolution of the source. 
- n_theta – The number of points in the theta direction. This value must be even. Default is 360. 
- n_check – For each image/projection find its common-lines with n_check images. If n_check is less than the total number of images, a random subset of n_check images is used. 
- max_shift – Determines maximum range for shifts for common-line detection as a proportion of the resolution. Default is 0.15. 
- shift_step – Resolution of shift estimation for common-line detection in pixels. Default is 1 pixel. 
- offsets_max_shift – Determines maximum range for shifts for 2D offset estimation as a proportion of the resolution. Default None inherits from max_shift. 
- offsets_shift_step – Resolution of shift estimation for 2D offset estimation in pixels. Default None inherits from shift_step. 
- hist_bin_width – Bin width in smoothing histogram (degrees). 
- full_width – Selection width around smoothed histogram peak (degrees). adaptive will attempt to automatically find the smallest number of `hist_bin_width`s required to find at least one valid image index. 
- mask – Option to mask src.images with a fuzzy mask (boolean). Default, True, applies a mask. 
 
 - estimate_rotations()¶
- Estimate rotation matrices using the common lines method with semi-definite programming. 
 
aspire.abinitio.commonline_sync module¶
- class aspire.abinitio.commonline_sync.CLSyncVoting(src, n_rad=None, n_theta=360, max_shift=0.15, shift_step=1, hist_bin_width=3, full_width=6, mask=True)¶
- Bases: - CLOrient3D,- SyncVotingMixin- Define a class to estimate 3D orientations using synchronization matrix and voting method. - The related publications are listed as below: Y. Shkolnisky, and A. Singer, Viewing Direction Estimation in Cryo-EM Using Synchronization, SIAM J. Imaging Sciences, 5, 1088-1110 (2012). - A. Singer, R. R. Coifman, F. J. Sigworth, D. W. Chester, Y. Shkolnisky, Detecting Consistent Common Lines in Cryo-EM by Voting, Journal of Structural Biology, 169, 312-322 (2010). - Initialize an object for estimating 3D orientations using synchronization matrix - Parameters:
- src – The source object of 2D denoised or class-averaged images with metadata 
- n_rad – The number of points in the radial direction 
- n_theta – The number of points in the theta direction. Default is 360. 
- max_shift – Determines maximum range for shifts as a proportion of the resolution. Default is 0.15. 
- shift_step – Resolution for shift estimation in pixels. Default is 1 pixel. 
- hist_bin_width – Bin width in smoothing histogram (degrees). 
- full_width – Selection width around smoothed histogram peak (degrees). adaptive will attempt to automatically find the smallest number of `hist_bin_width`s required to find at least one valid image index. 
- mask – Option to mask src.images with a fuzzy mask (boolean). Default, True, applies a mask. 
 
 - estimate_rotations()¶
- Estimate orientation matrices for all 2D images using synchronization matrix 
 - syncmatrix_vote()¶
- Construct the synchronization matrix using voting method - A pre-computed common line matrix is required as input. 
 
aspire.abinitio.commonline_sync3n module¶
- class aspire.abinitio.commonline_sync3n.CLSync3N(src, n_rad=None, n_theta=360, max_shift=0.15, shift_step=1, hist_bin_width=1, full_width='adaptive', epsilon=0.01, max_iters=1000, sigma=3, seed=None, mask=True, S_weighting=False, J_weighting=False, hist_intervals=100, disable_gpu=False)¶
- Bases: - CLOrient3D,- SyncVotingMixin- Define a class to estimate 3D orientations using common lines Sync3N methods (2017). - Ido Greenberg, Yoel Shkolnisky, Common lines modeling for reference free Ab-initio reconstruction in cryo-EM, Journal of Structural Biology, Volume 200, Issue 2, 2017, Pages 106-117, ISSN 1047-8477, https://doi.org/10.1016/j.jsb.2017.09.007. - Initialize object for estimating 3D orientations. - Parameters:
- src – The source object of 2D denoised or class-averaged images with metadata 
- n_rad – The number of points in the radial direction 
- n_theta – The number of points in the theta direction 
- max_shift – Maximum range for shifts as a proportion of box size. Default = 0.15. 
- shift_step – Step size of shift estimation in pixels. Default = 1 pixel. 
- hist_bin_width – Bin width in smoothing histogram (degrees). 
- full_width – Selection width around smoothed histogram peak (degrees). adaptive will attempt to automatically find the smallest number of `hist_bin_width`s required to find at least one valid image index. 
- epsilon – Tolerance for the power method. 
- max_iter – Maximum iterations for the power method. 
- sigma – Voting contribution smoothing factor. 
- seed – Optional seed for RNG. 
- mask – Option to mask src.images with a fuzzy mask (boolean). Default, True, applies a mask. 
- S_weighting – Optionally apply probabilistic weighting to the S matrix. 
- J_weighting – Optionally use J weights instead of signs when computing signs_times_v. 
- hist_intervals – Number of histogram bins used to compute triangle scores when S_weighting enabled. 
- disable_gpu – Disables GPU acceleration; forces CPU only code for this module. Defaults to automatically using GPU when available. 
 
 - estimate_rotations()¶
- Estimate rotation matrices. - Returns:
- Array of rotation matrices, size n_imgx3x3. 
 
 
aspire.abinitio.sync_voting module¶
- class aspire.abinitio.sync_voting.SyncVotingMixin¶
- Bases: - object- SyncVotingMixin is a mixin implementing methods for the synchronization voting algorithm which are shared by CLSynVoting and CLSymmetryC3C4