aspire.covariance package

Submodules

aspire.covariance.covar module

class aspire.covariance.covar.CovarianceEstimator(*args, **kwargs)

Bases: Estimator

An object representing a 2*L-by-2*L-by-2*L array containing the non-centered Fourier transform of the mean least-squares estimator kernel. Convolving a volume with this kernel is equal to projecting and backproject-ing that volume in each of the projection directions (with the appropriate amplitude multipliers and CTFs) and averaging over the whole dataset. Note that this is a non-centered Fourier transform, so the zero frequency is found at index 1.

Parameters:
  • srcImageSource to be used for estimation.

  • basis – 3D Basis to be used during estimation.

  • batch_size – Optional batch size of images drawn from src during back projection and kernel estimation steps.

  • preconditioner – Optional kernel preconditioner (string). Currently supported options are “circulant” or None.

  • checkpoint_iterations – Optionally save cg estimated Volume instance periodically each checkpoint_iterations. Setting to None disables, otherwise checks for positive integer.

  • checkpoint_prefix – Optional path prefix for cg checkpoint files. If the parent directory does not exist, creation is attempted. _iter{N} will be appended to the prefix.

  • maxiter – Optional max number of cg iterations before returning. This should be used in conjunction with checkpoint_iterations to prevent excessive disk usage. None disables.

  • boost – Option to use src symmetry to boost number of images used for mean estimation (Boolean). Default of True employs symmetry boosting.

apply_kernel(coef, kernel=None, packed=False)

Applies the kernel represented by convolution

Parameters:
  • coef – The volume matrix (6 dimensions) to be convolved (but see the packed argument below).

  • kernel – a Kernel object. If None, the kernel for this Estimator is used.

  • packed – whether the coef matrix represents an isometrically mapped packed vector, through the symmat_to_vec_iso function. In this case, the function expands coef into a symmetric matrix internally, and returns a packed vector in return.

Returns:

The result of evaluating coef in the given basis, convolving with the kernel given by kernel, and backprojecting into the basis. If packed is True, then the isometrically mapped packed vector is returned instead.

compute_kernel()
conj_grad(b_coef, tol=1e-05, regularizer=0)
estimate(mean_vol, noise_variance, tol=1e-05, regularizer=0)

Return an estimate as a Volume instance.

src_backward(mean_vol, noise_variance, shrink_method=None)

Apply adjoint mapping to source

Returns:

The sum of the outer products of the mean-subtracted images in src, corrected by the expected noise contribution and expressed as coefficients of basis.

aspire.covariance.covar2d module

class aspire.covariance.covar2d.BatchedRotCov2D(src, basis=None, batch_size=8192)

Bases: RotCov2D

Perform batchwise rotationally equivariant 2D covariance estimation from an ImageSource objects. This is done with a single pass through the data, processing moderately-sized batches one at a time. The rotational equivariance is achieved by decomposing images in a steerable Fourier–Bessel basis. For more information, see

T. Bhamre, T. Zhang, and A. Singer, “Denoising and covariance estimation of single particle cryo-EM images”, J. Struct. Biol. 195, 27-81 (2016). DOI: 10.1016/j.jsb.2016.04.013

Parameters:
  • src – The ImageSource object from which the sample images are to be extracted.

  • basis – The FBBasis2D object used to decompose the images. By default, this is set to FFBBasis2D((src.L, src.L)).

  • batch_size – The number of images to process at a time (default 8192).

constructor of an object for 2D covariance analysis

get_covar(noise_var=0, mean_coef=None, covar_est_opt=None, make_psd=True)

Calculate the block diagonal covariance matrix in the basis coefficients.

Parameters:
  • noise_var – The variance of the noise in the images (default 1)

  • mean_coef – If specified, overrides the mean coefficient vector used to calculate the covariance (default self.get_mean()).

  • :covar_est_opt

    The estimation parameters for obtaining the covariance matrix in the form of a dictionary. Keys include: - ‘shrinker’: The type of shrinkage we apply to the right-hand side

    in the normal equations. Can be ‘None’, in which case no shrinkage is performed. For a list of shrinkers, see the documentation of shrink_covar.

    • ’verbose’: Verbosity (integer) of the conjugate gradient algorithm (see documentation for conj_grad, default zero).

    • ’max_iter’: Maximum number of conjugate gradient iterations (see documentation for conj_grad, default 250).

    • ’iter_callback’: Callback performed at the end of an iteration (see documentation for conj_grad, default []).

    • ’store_iterates’: Determines whether to store intermediate iterates (see documentation for conj_grad, default False).

    • ’rel_tolerance’: Relative stopping tolerance of the conjugate gradient algorithm (see documentation for conj_grad, default 1e-12).

    • ’precision’: Precision of conjugate gradient algorithm (see documentation for conj_grad, default ‘float64’)

  • make_psd – If True, make the covariance matrix positive semidefinite

Returns:

The block diagonal matrix containing the basis coefficients (in self.basis) for the estimated covariance matrix. These are implemented using BlkDiagMatrix.

get_cwf_coefs(coefs, ctf_basis, ctf_idx, mean_coef, covar_coef, noise_var=0)

Estimate the expansion coefficients using the Covariance Wiener Filtering (CWF) method.

Parameters:
  • coefs – A coefficient vector (or an array of coefficient vectors) to be calculated.

  • ctf_basis – The CTF functions in the Basis expansion.

  • ctf_idx – An array of the CTF function indices for all 2D images. If ctf_basis or ctf_idx is None, the identity filter will be applied.

  • mean_coef – The mean value vector from all images.

  • covar_coef – The block diagonal covariance matrix of the clean coefficients represented by a cell array.

  • noise_var – The estimated variance of noise. The value should be zero for coefs from clean images of simulation data.

Returns:

The estimated coefficients of the unfiltered images in certain math basis. These are obtained using a Wiener filter with the specified covariance for the clean images and white noise of variance noise_var for the noise.

get_mean()

Calculate the rotationally invariant mean image in the basis coefficients.

Returns:

The mean coefficient vector in self.basis.

class aspire.covariance.covar2d.RotCov2D(basis)

Bases: object

Define a class for performing Cov2D analysis with CTF information described in

T. Bhamre, T. Zhang, and A. Singer, “Denoising and covariance estimation of single particle cryo-EM images”, J. Struct. Biol. 195, 27-81 (2016). DOI: 10.1016/j.jsb.2016.04.013

constructor of an object for 2D covariance analysis

get_covar(coefs, ctf_basis=None, ctf_idx=None, mean_coef=None, do_refl=True, noise_var=0, covar_est_opt=None, make_psd=True)

Calculate the covariance matrix from the expansion coefficients and CTF information.

Parameters:
  • coefs – A coefficient vector (or an array of coefficient vectors) to be calculated.

  • ctf_basis – The CTF functions in the Basis expansion.

  • ctf_idx – An array of the CTF function indices for all 2D images. If ctf_basis or ctf_idx is None, the identity filter will be applied.

  • mean_coef – The mean value vector from all images.

  • noise_var – The estimated variance of noise. The value should be zero for coefs from clean images of simulation data.

  • covar_est_opt – The optimization parameter list for obtaining the Cov2D matrix.

  • make_psd – If True, make the covariance matrix positive semidefinite

Returns:

The basis coefficients of the covariance matrix in the form of cell array representing a block diagonal matrix. These block diagonal matrices are implemented as BlkDiagMatrix instances. The covariance is calculated from the images represented by the coefs array, along with all possible rotations and reflections. As a result, the computed covariance matrix is invariant to both reflection and rotation. The effect of the filters in ctf_basis are accounted for and inverted to yield a covariance estimate of the unfiltered images.

get_cwf_coefs(coefs, ctf_basis=None, ctf_idx=None, mean_coef=None, covar_coef=None, noise_var=0)

Estimate the expansion coefficients using the Covariance Wiener Filtering (CWF) method.

Parameters:
  • coefs – A coefficient vector (or an array of coefficient vectors) to be calculated.

  • ctf_basis – The CTF functions in the Basis expansion.

  • ctf_idx – An array of the CTF function indices for all 2D images. If ctf_basis or ctf_idx is None, the identity filter will be applied.

  • mean_coef – The mean value vector from all images.

  • covar_coef – The block diagonal covariance matrix of the clean coefficients represented by a cell array.

  • noise_var – The estimated variance of noise. The value should be zero for coefs from clean images of simulation data.

Returns:

The estimated coefficients of the unfiltered images in certain math basis. These are obtained using a Wiener filter with the specified covariance for the clean images and white noise of variance noise_var for the noise.

get_mean(coefs, ctf_basis=None, ctf_idx=None)

Calculate the mean vector from the expansion coefficients with CTF information.

Parameters:
  • coefs – A coefficient vector (or an array of coefficient vectors) to be averaged.

  • ctf_basis – The CTF functions in the Basis expansion.

  • ctf_idx – An array of the CTF function indices for all 2D images. If ctf_basis or ctf_idx is None, the identity filter will be applied.

Returns:

The mean value vector for all images.

shrink_covar_backward(b, b_noise, n, noise_var, shrinker)

Apply the shrinking method to the 2D covariance of coefficients.

Parameters:
  • b – An input coefficient covariance.

  • b_noise – The noise covariance.

  • noise_var – The estimated variance of noise.

  • shrinker – The shrinking method.

Returns:

The shrinked 2D covariance coefficients.

aspire.covariance.covar2d.shrink_covar(covar, noise_var, gamma, shrinker='frobenius_norm')

Shrink the covariance matrix

Parameters:
  • covar_in – An input covariance matrix

  • noise_var – The estimated variance of noise

  • gamma – An input parameter to specify the maximum values of eigen values to be neglected.

  • shrinker – An input parameter to select different shrinking methods.

Returns:

The shrinked covariance matrix

Module contents