aspire.basis package

Submodules

aspire.basis.basis module

class aspire.basis.basis.Basis(size, ell_max=None, dtype=<class 'numpy.float32'>)

Bases: object

Define a base class for expanding 2D particle images and 3D structure volumes

Initialize an object for the base of basis class

Parameters:
  • size – The size of the vectors for which to define the basis. Currently only square images and cubic volumes are supported.

  • ell_max – The maximum order ell of the basis elements. If no input (= None), it will be set to np.Inf and the basis includes all ell such that the resulting basis vectors are concentrated below the Nyquist frequency (default Inf).

evaluate(v)

Evaluate coefficient vector in basis

Parameters:

vCoef instance containing the coefficients to be evaluated. The first dimension must correspond to the number of coefficient vectors, while the second must correspond to self.count.

Returns:

The evaluation of the coefficient vector(s) v for this basis. This is an Image or a Volume object containing one image/volume for each coefficient vector, and of size self.sz.

evaluate_t(v)

Evaluate coefficient in dual basis

Parameters:

v – An Image or Volume object whose size matches self.sz.

Returns:

The evaluation of the Image or Volume object v in the dual basis of basis. This is an array of vectors whose first dimension equals the number of images/volumes in v. and whose second dimension is self.count.

expand(x, tol=None, atol=0)

Obtain coefficients in the basis from those in standard coordinate basis

This is a similar function to evaluate_t but with more accuracy by using the cg optimizing of linear equation, Ax=b.

Parameters:
  • x – An array whose last two or three dimensions are to be expanded the desired basis. These dimensions must equal self.sz.

  • tol – Relative tolerance for convergence, norm(residual) <= max(tol*norm(b), atol). Deafult None sets to dtype’s `eps`*10.

  • atol – Absolute tolerance for convergence, norm(residual) <= max(tol*norm(b), atol).

Returns:

The coefficients of v expanded in the desired basis. The last dimension of v is with size of count and the first dimensions of the return value correspond to those first dimensions of x.

mat_evaluate(V)

Evaluate coefficient matrix in basis

Parameters:

V – A coefficient matrix of size self.count-by- self.count to be evaluated.

Returns:

A multidimensional matrix of size self.sz-by -self.sz corresponding to the evaluation of V in this basis.

mat_evaluate_t(X)

Evaluate coefficient matrix in dual basis

Parameters:

X – The coefficient array of size self.sz-by-self.sz to be evaluated.

Returns:

The evaluation of X in the dual basis. This is self.count-by-self.count. matrix. If V is a matrix of size self.count-by-self.count, B is the change-of-basis matrix of basis, and x is a multidimensional matrix of size basis.sz-by-basis.sz, the function calculates V = B’ * X * B, where the rows of B, rows of ‘X’, and columns of X are read as vectorized arrays.

norms()

Calculate the normalized factors of basis functions

class aspire.basis.basis.Coef(basis, data, dtype=None)

Bases: object

Numpy interoperable container for stacks of real coefficient vectors. Each Coef instance has an associated Basis.

A stack of one or more coefficient arrays.

The stack can be multidimensional with stack_size equal to the product of the stack dimensions. Singletons will be expanded into a 1D stack of length one.

The last axes always represents the coefficient count.

Parameters:
  • basisBasis associated with data coefficients.

  • data – Numpy array containing image data with shape (…, count).

  • dtype – Optionally cast data to this dtype. Defaults to data.dtype.

Returns:

Coef instance holding data.

asnumpy()

Return image data as a (<stack>, count) read-only array view.

Returns:

read-only ndarray view

by_indices(**kwargs)

Select coefficients by indices (radial, angular).

See SteerableBasis.indices_mask for argument details.

Returns:

Numpy array.

copy()

Return a new Coef instance with a deep copy of the data.

evaluate()

Return the evaluation of coefficients in the associated basis.

rotate(radians, refl=None)

Returns coefs rotated counter-clockwise by radians.

Raises error if underlying coef basis does not support rotations.

Parameters:
  • radians – Rotation in radians.

  • refl – Optional reflect image (about y=0) (bool)

Returns:

rotated coefs.

shift(shifts)

Returns coefs shifted by shifts.

This will transform to real cartesian space, shift, and transform back to basis space.

Parameters:
  • coef – Basis coefs.

  • shifts – Shifts in pixels (x,y). Shape (1,2) or (len(coef), 2).

Returns:

coefs of shifted images.

property size

Return np.size of underlying data.

This should be stack_size * count, or len(self) * count.

stack_reshape(*args)

Reshape the stack axis.

*argsargs:

Integer(s) or tuple describing the intended shape.

Returns:

Coef instance

to_complex()

Convert and return real coefficients as ComplexCoef.

to_real()

Not implemented for real Coef.

class aspire.basis.basis.ComplexCoef(basis, data, dtype=None)

Bases: Coef

Numpy interoperable container for stacks of complex coefficient vectors. Each ComplexCoef instance has an associated Basis.

A stack of one or more coefficient arrays.

The stack can be multidimensional with stack_size equal to the product of the stack dimensions. Singletons will be expanded into a 1D stack of length one.

The last axes always represents the coefficient count.

Parameters:
  • basisBasis associated with data coefficients.

  • data – Numpy array containing image data with shape (…, count).

  • dtype – Optionally cast data to this dtype. Defaults to data.dtype.

Returns:

Coef instance holding data.

evaluate()

Return the evaluation of coefficients in the associated basis.

rotate(radians, refl=None)

Returns coefs rotated counter-clockwise by radians.

Raises error if underlying coef basis does not support rotations.

Parameters:
  • radians – Rotation in radians.

  • refl – Optional reflect image (about y=0) (bool)

Returns:

Rotated ComplexCoefs.

shift(shifts)

Returns complex coefs shifted by shifts.

This will transform to real cartesian space, shift, and transform back to basis space.

Parameters:
  • coef – Basis coefs.

  • shifts – Shifts in pixels (x,y). Shape (1,2) or (len(coef), 2).

Returns:

Complex coefs of shifted images.

to_complex()

Not implemented for ComplexCoef.

to_real()

Convert and return complex coefficients as Coef.

aspire.basis.basis_utils module

Define related utility functions for Fourier–Bessel (2D), Spherical Fourier–Bessel (3D) and prolate spheroidal wave function (PSWF) objects.

aspire.basis.basis_utils.all_besselj_zeros(ell, r)

Compute the zeros of the order ell Bessel function which are less than r.

Parameters:
  • ell – The real number order of the Bessel function.

  • r – The upper bound for zeros returned.

Return n, r0:

The number of zeros and the zeros themselves as a NumPy array.

aspire.basis.basis_utils.besselj_newton(nu, z0, max_iter=10)

Uses the Newton-Raphson method to compute the zero(s) of the Bessel function with order nu with initial guess(es) z0.

Parameters:
  • nu – The real number order of the Bessel function.

  • z0 – (Array-like) The initial guess(es) for the root-finding algorithm.

  • max_iter – Maximum number of iterations for Newton-Raphson (default: 10).

Return z:

(Array-like) The estimated root(s).

aspire.basis.basis_utils.besselj_zeros(nu, k)

Finds the first k zeros of the Bessel function of order nu, i.e. J_nu. Adapted from “zerobess.m” by Jonas Lundgren <splinefit@gmail.com>

Parameters:
  • nu – The real number order of the Bessel function (must be positive and <1e7).

  • k – The number of zeros to return (must be >= 3).

Return z:

A 1D NumPy array of the first k zeros.

aspire.basis.basis_utils.check_besselj_zeros(nu, z)

Sanity-check a sequence of estimated zeros of the Bessel function with order nu.

Parameters:
  • nu – The real number order of the Bessel function.

  • z – (Array-like) A sequence of postulated zeros.

Return result:

True or False.

aspire.basis.basis_utils.d_decay_approx_fun(a, b, c, d)
aspire.basis.basis_utils.k_operator(nu, x)
aspire.basis.basis_utils.lgwt(ndeg, a, b, dtype=<class 'numpy.float32'>)

Compute Legendre-Gauss quadrature

Generates the Legendre-Gauss nodes and weights on an interval [a, b] with truncation order of ndeg for computing definite integrals using Legendre-Gauss quadrature. Suppose you have a continuous function f(x) which is defined on [a, b] which you can evaluate at any x in [a, b]. Simply evaluate it at all of the values contained in the x vector to obtain a vector f, then compute the definite integral using sum(f.*w);

This is a 2rapper for numpy.polynomial leggauss which outputs only in the range of (-1, 1).

Parameters:
  • ndeg – truncation order, that is, the number of nodes.

  • b (a,) – The endpoints of the interval over which the quadrature is defined.

Return x, w:

The quadrature nodes and weights.

aspire.basis.basis_utils.norm_assoc_legendre(j, m, x)

Evaluate the normalized associated Legendre polynomial

Parameters:
  • j – The order of the associated Legendre polynomial, must satisfy |m| < j.

  • m – The degree of the associated Legendre polynomial, must satisfy |m| < j.

  • x – An array of values between -1 and +1 on which to evaluate.

Returns:

The normalized associated Legendre polynomial evaluated at corresponding x.

aspire.basis.basis_utils.p_n(n, alpha, beta, x)

The first n jacobi polynomial of x as defined in Yoel’s PSWF paper, eq (2), page 6

Parameters:
  • n – int, > 0 Number of polynomials to compute

  • alpha – float, > -1

  • beta – float, > -1

  • x – (m,) ndarray

Returns:

v: (m, n + 1) ndarray v[:, i] = P^{(alpha, beta)}_n(x) as defined in the paper

aspire.basis.basis_utils.real_sph_harmonic(j, m, theta, phi)

Evaluate a real spherical harmonic

Parameters:
  • j – The order of the spherical harmonic. These must satisfy |m| < j.

  • m – The degree of the spherical harmonic. These must satisfy |m| < j.

  • theta – The spherical coordinates of the points at which we want to evaluate the real spherical harmonic. theta is the latitude between 0 and pi

  • phi – The spherical coordinates of the points at which we want to evaluate the real spherical harmonic. phi is the longitude, between 0 and 2*pi

Returns:

The real spherical harmonics evaluated at the points (theta, phi).

aspire.basis.basis_utils.sph_bessel(ell, r)

Compute spherical Bessel function values.

Parameters:
  • ell – The order of the spherical Bessel function.

  • r – The coordinates where the function is to be evaluated.

Returns:

The value of j_ell at r.

aspire.basis.basis_utils.t_radial_part_mat(x, n, j, m)
aspire.basis.basis_utils.t_x_derivative_mat(t1, t2, x, big_n, range_array, approx_length)
aspire.basis.basis_utils.t_x_mat(x, n, j, approx_length)
aspire.basis.basis_utils.t_x_mat_dot(x, n, j, approx_length)
aspire.basis.basis_utils.unique_coords_nd(N, ndim, shifted=False, normalized=True, dtype=<class 'numpy.float32'>)

Generate unique polar coordinates from 2D or 3D rectangular coordinates.

Parameters:
  • N – length size of a square or cube.

  • ndim – number of dimension, 2 or 3.

  • shifted – shifted half pixel or not for odd N.

  • normalized – normalize the grid or not.

Returns:

The unique polar coordinates in 2D or 3D

aspire.basis.fb module

class aspire.basis.fb.FBBasisMixin

Bases: object

FBBasisMixin is a mixin implementing methods specific to Fourier-Bessel expansions, to be inherited by Fourier-Bessel subclasses of Basis.

aspire.basis.fb_2d module

class aspire.basis.fb_2d.FBBasis2D(size, ell_max=None, dtype=<class 'numpy.float32'>)

Bases: SteerableBasis2D, FBBasisMixin

Define a derived class using the Fourier-Bessel basis for mapping 2D images

The expansion coefficients of 2D images on this basis are obtained by the least squares method. The algorithm is described in the publication: Z. Zhao, A. Singer, Fourier-Bessel Rotational Invariant Eigenimages, The Journal of the Optical Society of America A, 30 (5), pp. 871-877 (2013).

Initialize an object for the 2D Fourier-Bessel basis class

Parameters:

size – The size of the vectors for which to define the basis. May be a 2-tuple or an integer, in which case a square basis is assumed. Currently only square images are supported.

Ell_max:

The maximum order ell of the basis elements. If no input (= None), it will be set to np.Inf and the basis includes all ell such that the resulting basis vectors are concentrated below the Nyquist frequency (default Inf).

basis_norm_2d(ell, k)

Calculate the normalized factors from radial and angular parts of a specified basis function

calculate_bispectrum(coef, flatten=False, filter_nonzero_freqs=False, freq_cutoff=None)

Calculate bispectrum for a set of coefs in this basis.

The Bispectum matrix is of shape:

(count, count, unique_radial_indices)

where count is the number of complex coefficients.

Parameters:
  • coef – Coefficients representing a (single) image expanded in this basis.

  • flatten – Optionally extract symmetric values (tril) and then flatten.

  • freq_cutoff – Truncate (zero) high k frequecies above (int) value, defaults off (None).

Returns:

Bispectum matrix (complex valued).

filter_to_basis_mat(*args, **kwargs)

See SteerableBasis2D.filter_to_basis_mat.

norms()

Calculate the normalized factors of basis functions

aspire.basis.fb_3d module

class aspire.basis.fb_3d.FBBasis3D(size, ell_max=None, dtype=<class 'numpy.float32'>)

Bases: Basis, FBBasisMixin

Define a derived class for direct spherical Harmonics Bessel basis expanding 3D volumes

# TODO: Methods that return dictionaries should return useful objects instead

Initialize an object for the 3D Fourier-Bessel basis class

Parameters:
  • size – The size of the vectors for which to define the basis. May be a 3-tuple or an integer, in which case a cubic basis is assumed. Currently only cubic images are supported.

  • ell_max – The maximum order ell of the basis elements. If no input (= None), it will be set to np.Inf and the basis includes all ell such that the resulting basis vectors are concentrated below the Nyquist frequency (default Inf).

basis_norm_3d(ell, k)

Calculate the normalized factor of a specified basis function.

indices()

Create the indices for each basis function

norms()

Calculate the normalized factors of basis functions

aspire.basis.ffb_2d module

class aspire.basis.ffb_2d.FFBBasis2D(size, ell_max=None, dtype=<class 'numpy.float32'>)

Bases: FBBasis2D

Define a derived class for Fast Fourier Bessel expansion for 2D images

The expansion coefficients of 2D images on this basis are obtained by a fast method instead of the least squares method. The algorithm is described in the publication: Z. Zhao, Y. Shkolnisky, A. Singer, Fast Steerable Principal Component Analysis, IEEE Transactions on Computational Imaging, 2 (1), pp. 1-12 (2016).​

Initialize an object for the 2D Fourier-Bessel basis class

Parameters:

size – The size of the vectors for which to define the basis. May be a 2-tuple or an integer, in which case a square basis is assumed. Currently only square images are supported.

Ell_max:

The maximum order ell of the basis elements. If no input (= None), it will be set to np.Inf and the basis includes all ell such that the resulting basis vectors are concentrated below the Nyquist frequency (default Inf).

filter_to_basis_mat(f, **kwargs)

See SteerableBasis2D.filter_to_basis_mat.

aspire.basis.ffb_3d module

class aspire.basis.ffb_3d.FFBBasis3D(size, ell_max=None, dtype=<class 'numpy.float32'>)

Bases: FBBasis3D

Define a derived class for fast spherical Harmonics Bessel basis expanding 3D volumes

# TODO: Methods that return dictionaries should return useful objects instead

Initialize an object for the 3D Fourier-Bessel basis class

Parameters:
  • size – The size of the vectors for which to define the basis. May be a 3-tuple or an integer, in which case a cubic basis is assumed. Currently only cubic images are supported.

  • ell_max – The maximum order ell of the basis elements. If no input (= None), it will be set to np.Inf and the basis includes all ell such that the resulting basis vectors are concentrated below the Nyquist frequency (default Inf).

aspire.basis.fle_2d module

class aspire.basis.fle_2d.FLEBasis2D(size, bandlimit=None, epsilon=1e-10, dtype=<class 'numpy.float32'>, match_fb=True)

Bases: SteerableBasis2D, FBBasisMixin

Define a derived class for Fast Fourier Bessel 2D expansion using interpolation

from Chebyshev nodes.

The algorithms used are described in the following publication: N. F. Marshall, O. Mickelin, A. Singer, Fast Expansion into Harmonics on the Disk:

A Steerable Basis with Fast Radial Convolution. (submitted)

https://arxiv.org/pdf/2207.13674.pdf

Parameters:
  • size – The size of the vectors for which to define the FLE basis. Currently only square images are supported.

  • bandlimit – Maximum frequency band for computing basis functions. Defaults to the resolution of the basis.

  • epsilon – Relative precision between FLE fast method and dense matrix multiplication.

  • dtype – Datatype of images and coefficients represented.

  • match_fb – This flag constructs basis functions identical to FBBasis2D. The initial heuristic for the number of basis functions, based on the image size, will be set to that of FBBasis2D, and the FLE frequency thresholding procedure to reduce the number of functions will not be carried out. This means the number of basis functions for a given image size will be identical across the two bases.

filter_to_basis_mat(f, **kwargs)

See SteerableBasis2D.filter_to_basis_mat.

lowpass(coefs, bandlimit)

Apply a low-pass filter to FLE coefficients coefs with threshold bandlimit.

Parameters:
  • coefs – A Coef instance containing FLE coefficients.

  • bandlimit – Integer bandlimit (max frequency).

Returns:

Band-limited coefficient array.

matrix_type

alias of DiagMatrix

radial_convolve(coefs, radial_img)

Convolve a stack of FLE coefficients with a 2D radial function.

Parameters:
  • coefs – A Coef instance containing FLE coefficients.

  • radial_img – A 2D NumPy array of size (self.nres, self.nres).

Returns:

Convolved FLE coefficients.

aspire.basis.fle_2d_utils module

aspire.basis.fle_2d_utils.barycentric_interp_sparse(target_points, known_points, numsparse)
Returns the sparse matrices that perform barycentric interpolation to compute values

of Betas at the points target_points at known points known_points, and the transpose of this operation. For each target point in target_points, only numsparse centered source points from known_points around the target point are used.

Performed via the method described in

“Barycentric Lagrange Interpolation”, Jean-Paul Berrut and Lloyd Trefethen. SIAM Review 2004 46:3, 501-517 https://people.maths.ox.ac.uk/trefethen/barycentric.pdf

Parameters:
  • target_points – The target set of points at which to evaluate the functions.

  • known_points – The points at which the values of the functions are known.

  • numsparse – Number of points used for interpolation around each target point.

Returns:

The interpolation matrix and its transpose as a 2-tuple.

aspire.basis.fle_2d_utils.precomp_transform_complex_to_real(ells)
Returns a sparse matrix that transforms coefficients into the complex

representation of the basis to coefficients in the real representation of the basis. See Remark 1.1 of Marshall, Mickelin, and Singer.

Parameters:

ells – The list of integer Bessel function orders.

Returns:

Sparse complex to real transformation matrix.

aspire.basis.fle_2d_utils.transform_complex_to_real(B, ells)
Transforms coefficients of the matrix B (see Eq. 3) from complex

to real. B is the linear transformation that takes FLE coefficients to images.

Parameters:
  • B – Complex matrix B.

  • ells – List of ells (Bessel function orders) in this basis.

Returns:

Transformed matrix.

aspire.basis.fpswf_2d module

class aspire.basis.fpswf_2d.FPSWFBasis2D(size, gamma_truncation=1.0, beta=1.0, dtype=<class 'numpy.float32'>)

Bases: PSWFBasis2D

Define a derived class for fast Prolate Spheroidal Wave Function (PSWF) expanding 2D images

The numerical evaluation for 2D PSWFs at arbitrary points in the unit disk is based on the fast method described in the papers as below:

1) Boris Landa and Yoel Shkolnisky, “Steerable principal components for space-frequency localized images”, SIAM J. Imag. Sci. 10, 508-534 (2017). 2) Boris Landa and Yoel Shkolnisky, “Approximation scheme for essentially bandlimited and space-concentrated functions on a disk”, Appl. Comput. Harmon. Anal. 43, 381-403 (2017). 3) Yoel Shkolnisky, “Prolate spheroidal wave functions on a disc-Integration and approximation of two-dimensional bandlimited functions”, Appl. Comput. Harmon. Anal. 22, 235-256 (2007).

Initialize an object for 2D prolate spheroidal wave function (PSWF) basis expansion using fast method.

Parameters:
  • size – The size of the vectors for which to define the basis and the image resolution. May be a 2-tuple or an integer, in which case a square basis is assumed. Currently only square images are supported.

  • gamma_trunc – Truncation parameter of PSWFs, between 0 and 1e6, which controls the length of the expansion and the approximation error. Smaller values (close to zero) guarantee smaller errors, yet longer expansions, and vice-versa. Note: Due to numerical considerations, do not exceed 1e6.

  • beta – Bandlimit ratio relative to the Nyquist rate, between 0 and 1. In general, the bandlimit is c = beta*pi*(size[0]//2), therefore for the default value beta = 1 there is no oversampling assumed. This parameter controls the bandlimit of the PSWFs.

  • dtype – Internal ndarray datatype.

filter_to_basis_mat(*args, **kwargs)

See SteerableBasis2D.filter_to_basis_mat.

aspire.basis.fpswf_3d module

class aspire.basis.fpswf_3d.FPSWFBasis3D(size, ell_max=None, dtype=<class 'numpy.float32'>)

Bases: PSWFBasis3D

Initialize an object for the base of basis class

Parameters:
  • size – The size of the vectors for which to define the basis. Currently only square images and cubic volumes are supported.

  • ell_max – The maximum order ell of the basis elements. If no input (= None), it will be set to np.Inf and the basis includes all ell such that the resulting basis vectors are concentrated below the Nyquist frequency (default Inf).

aspire.basis.fspca module

class aspire.basis.fspca.FSPCABasis(src, basis=None, noise_var=None, components=None, batch_size=512)

Bases: SteerableBasis2D

A class for Fast Steerable Principal Component Analaysis basis.

FSPCA is an extension to Fourier Bessel representations (provided asF BBasis2D/FFBBasis2D), which computes combinations of basis coefficients coresponding to the princicpal components of image(s) represented in the provided basis.

The principal components are computed from eigen decomposition of the covariance matrix, and when evaluated into the real domain and reshaped form the set of eigenimages.

The algorithm is described in the publication: Z. Zhao, Y. Shkolnisky, A. Singer, Fast Steerable Principal Component Analysis, IEEE Transactions on Computational Imaging, 2 (1), pp. 1-12 (2016).

Parameters:
  • src – Source instance

  • basis – Optional Fourier Bessel Basis (usually FFBBasis2D)

  • components – Optionally assign number of principal components to use for the FSPCA basis. Default value of None will use self.basis.count.

  • noise_var – Optionally assign noise variance. Default value of None will estimate noise with WhiteNoiseEstimator. Use 0 when using clean images so cov2d skips applying noisy covar coefs..

  • batch_size – Batch size for computing basis coefficients. batch_size is also passed to BatchedRotCov2D.

build()

Computes the FSPCA basis.

This may take some time for large image stacks.

calculate_bispectrum(coef, flatten=False, filter_nonzero_freqs=False, freq_cutoff=None)

Calculate bispectrum for a set of coefs in this basis.

The Bispectum matrix is of shape:

(count, count, unique_radial_indices)

where count is the number of complex coefficients.

Parameters:
  • coef – Coefficients representing a (single) image expanded in this basis.

  • flatten – Optionally extract symmetric values (tril) and then flatten.

  • filter_nonzero_freqs – Remove indices corresponding to zero frequency (defaults False).

  • freq_cutoff – Truncate (zero) high k frequecies above (int) value, defaults off (None).

Returns:

Bispectum matrix (complex valued).

eigen_images()

Return the eigen images of the FSPCA basis, evaluated to image space.

This may be used to implot visualizations of the eigenvectors.

Ordering corresponds to FSPCA eigvals.

property eigvals

Return the eigenvals of FSPCABasis as Numpy array.

evaluate(c)

Take FSPCA coefs and evaluate to Fourier Bessel (self.basis) ceofs.

Parameters:

c – Stack of coefs in the FSPCABasis to be evaluated.

Returns:

The (real) coefs representing a stack of images in self.basis

evaluate_to_image_basis(c)

Take FSPCA coefs and evaluate as image in the standard coordinate basis.

Parameters:

c – Stack of coefs in the FSPCABasis to be evaluated.

Returns:

The Image instance representing a stack of images in the standard 2D coordinate basis..

expand(x)

Take a Fourier-Bessel coefs and express as FSPCA coefs.

Note each FSPCA coef corresponds to a linear combination Fourier Bessel basis vectors, described by an eigenvector in FSPCA.

Parameters:

x – Coefs representing a stack in the Fourier Bessel basis.

Returns:

Stack of coefs in the FSPCABasis.

expand_from_image_basis(x)

Take an image in the standard coordinate basis and express as FSPCA coefs.

Note each FSPCA coef corresponds to a linear combination Fourier Bessel basis vectors, described by an eigenvector in FSPCA.

Parameters:

x – The Image instance representing a stack of images in the standard 2D coordinate basis to be evaluated.

Returns:

Stack of coefs in the FSPCABasis.

filter_to_basis_mat(f)

Convert a filter into a basis representation.

Parameters:

fFilter object, usually a CTFFilter.

Returns:

Representation of filter in basis. Return type will be based on the class’s matrix_type.

shift(coef, shifts)

Returns coefs shifted by shifts.

This will transform to real cartesian space, shift, and transform back to Polar Fourier-Bessel space.

Parameters:
  • coef – Basis coefs.

  • shifts – Shifts in pixels (x,y). Shape (1,2) or (len(coef), 2).

Returns:

coefs of shifted images.

to_complex(coef)

Return complex valued representation of coefficients. This can be useful when comparing or implementing methods from literature.

There is a corresponding method, to_real.

Parameters:

coef – Coefficients from this basis.

Returns:

Complex coeficent representation from this basis.

to_real(complex_coef)

Return real valued representation of complex coefficients. This can be useful when comparing or implementing methods from literature.

There is a corresponding method, to_complex.

Parameters:

complex_coef – Complex coefficients from this basis.

Returns:

Real coefficient representation from this basis.

aspire.basis.pswf_2d module

class aspire.basis.pswf_2d.PSWFBasis2D(size, gamma_trunc=1.0, beta=1.0, dtype=<class 'numpy.float32'>)

Bases: SteerableBasis2D

Define a derived class for direct Prolate Spheroidal Wave Function (PSWF) expanding 2D images

The numerical evaluation for 2D PSWFs at arbitrary points in the unit disk is based on the direct method described in the papers as below:

1) Boris Landa and Yoel Shkolnisky, “Steerable principal components for space-frequency localized images”, SIAM J. Imag. Sci. 10, 508-534 (2017). 2) Boris Landa and Yoel Shkolnisky, “Approximation scheme for essentially bandlimited and space-concentrated functions on a disk”, Appl. Comput. Harmon. Anal. 43, 381-403 (2017). 3) Yoel Shkolnisky, “Prolate spheroidal wave functions on a disc-Integration and approximation of two-dimensional bandlimited functions”, Appl. Comput. Harmon. Anal. 22, 235-256 (2007).

Initialize an object for 2D PSWF basis expansion using direct method

Parameters:
  • size – The size of the vectors for which to define the basis and the image resolution. May be a 2-tuple or an integer, in which case a square basis is assumed. Currently only square images are supported.

  • gamma_trunc – Truncation parameter of PSWFs, between 0 and 1e6, which controls the length of the expansion and the approximation error. Smaller values (close to zero) guarantee smaller errors, yet longer expansions, and vice-versa. Note: Due to numerical considerations, do not exceed 1e6.

  • beta – Bandlimit ratio relative to the Nyquist rate, between 0 and 1. In general, the bandlimit is c = beta*pi*(size[0]//2), therefore for the default value beta = 1 there is no oversampling assumed. This parameter controls the bandlimit of the PSWFs.

  • dtype – Internal ndarray datatype.

filter_to_basis_mat(*args, **kwargs)

See SteerableBasis2D.filter_to_basis_mat.

matrix_type

alias of BlkDiagMatrix

pswf_func2d(big_n, n, bandlimit, phi_approximate_error, r, w)

Calculate the eigenvalues and eigenvectors of PSWF basis functions for all N’s and n’s

Parameters:
  • big_n – The integer N in PSWF basis.

  • n – The integer n in PSWF basis.

  • bandlimit – The band limit estimated by beta * pi * rcut.

  • phi_approximate_error – The input approximate error for phi.

  • r – The Legendre–Gauss quadrature nodes.

  • w – The Legendre–Gauss quadrature weights.

Returns:

alpha_n (ndarray): the eigen-values for N. d_vec (ndarray): the corresponding eigen-vectors for alpha_n. approx_length (int): the number of eigenvalues,len(alpha_n).

aspire.basis.pswf_3d module

class aspire.basis.pswf_3d.PSWFBasis3D(size, ell_max=None, dtype=<class 'numpy.float32'>)

Bases: Basis

Initialize an object for the base of basis class

Parameters:
  • size – The size of the vectors for which to define the basis. Currently only square images and cubic volumes are supported.

  • ell_max – The maximum order ell of the basis elements. If no input (= None), it will be set to np.Inf and the basis includes all ell such that the resulting basis vectors are concentrated below the Nyquist frequency (default Inf).

aspire.basis.pswf_utils module

class aspire.basis.pswf_utils.BNMatrix(big_n, band_limit, approx_length, dtype=<class 'numpy.float32'>)

Bases: object

Define a class to compute the B_N matrix ( with elements of b^N_mn) denoting the matrix of the operator L_c with respect to the basis T_Nn(x) by elements of b^N_mn. The matrix element b^N_mn is calculated as < T_Nm, L_c T_Nn > as shown in the paper below:

Yoel Shkolnisky, “Prolate spheroidal wave functions on a disc-Integration and approximation of two-dimensional bandlimited functions”, Appl. Comput. Harmon. Anal. 22, 235-256 (2007).

Initial an object to compute the B_N matrix ( with elements of b^N_mn).

Parameters:
  • big_n – A positive integer represented by N.

  • band_limit – The band limit.

  • approx_length – The approximated length represented by n.

dense_mat()

Represent the full B_N matrix by a 2D array.

Returns:

mat: (M,M), ndarray The full BN matrix.

get_eig_vectors()

Calculated the eigenvalues and eigenvectors of B_N matrix.

Returns:

v: (M,M) ndarray

The normalized eigenvectors corresponding to the eigenvalues, v[:, i] is corresponding to the w[i]. In each eigenvector v[:, i], v[argmax(abs(v[:, i])), i] >= 0.

w: (M,) ndarray

The eigenvalues in descending order.

shape()

Ruturn the shape of B_N matrix.

Returns:

tuple: (n, n) The dense matrix shape

aspire.basis.steerable module

class aspire.basis.steerable.SteerableBasis2D(*args, **kwargs)

Bases: Basis, ABC

SteerableBasis2D is an extension of Basis that is expected to have rotation (steerable) and calculate_bispectrum methods.

Initialize an object for the base of basis class

Parameters:
  • size – The size of the vectors for which to define the basis. Currently only square images and cubic volumes are supported.

  • ell_max – The maximum order ell of the basis elements. If no input (= None), it will be set to np.Inf and the basis includes all ell such that the resulting basis vectors are concentrated below the Nyquist frequency (default Inf).

property blk_diag_cov_shape

Return the BlkDiagMatrix partition shapes.

If the shape has already been cached, returns cached value. Otherwise, will compute the shape and cache in this instance.

calculate_bispectrum(complex_coef, flatten=False, filter_nonzero_freqs=False, freq_cutoff=None)

Calculate bispectrum for a set of coefs in this basis.

The Bispectum matrix is of shape:

(count, count, unique_radial_indices)

where count is the number of complex coefficients.

Parameters:
  • coef – Coefficients representing a (single) image expanded in this basis.

  • flatten – Optionally extract symmetric values (tril) and then flatten.

  • filter_nonzero_freqs – Remove indices corresponding to zero frequency (defaults False).

  • freq_cutoff – Truncate (zero) high k frequecies above (int) value, defaults off (None).

Returns:

Bispectum matrix (complex valued).

complex_rotate(complex_coef, radians, refl=None)

Returns complex coefs rotated counter-clockwise by radians.

This implementation uses the complex exponential. It is kept in the code for documentation and reference purposes.

To invoke in code:

self.to_real(
self.complex_rotate(

self.to_complex(coef), radians, refl)

)

)

Parameters:
  • complex_coef – Basis coefs (in complex representation).

  • radians – Rotation in radians.

  • refl – Optional reflect image (about y=0) (bool)

Returns:

rotated (complex) coefs.

abstract filter_to_basis_mat(f, method='evaluate_t', truncate=True)

Convert a filter into a basis operator representation.

Parameters:
  • fFilter object, usually a CTFFilter.

  • methodevaluate_t or expand.

  • truncate – Optionally, truncate dense matrix to BlkDiagMatrix. Defaults to True.

Returns:

Representation of filter as basis operator. Return type will be based on the class’s matrix_type.

indices_mask(**kwargs)

Given radial= or angular= expressions, return (count,) shaped mask where values satisfying the expression are True.

Examples

No args yield all indices. angular=0 creates a mask for selecting coefficients with zero angular indices. `angular=1, radial=2 selects coefficients satisfying angular index of 1 _and_ radial index of 2. More advanced operations can combine indices attributes.

angular=self.angular_indices>=0, radial=r selects coefficients with non negative angular indices and some radial index r.

Returns:

Boolen mask of shape (count,). Intended to be broadcast with Coef containers.

matrix_type

alias of BlkDiagMatrix

rotate(coef, radians, refl=None)

Returns coefs rotated counter-clockwise by radians.

Parameters:
  • coef – Basis coefs.

  • radians – Rotation in radians.

  • refl – Optional reflect image (about y=0) (bool)

Returns:

rotated coefs.

shift(coef, shifts)

Returns coefs shifted by shifts.

This will transform to real cartesian space, shift, and transform back to Polar Fourier-Bessel space.

Parameters:
  • coef – Basis coefs.

  • shifts – Shifts in pixels (x,y). Shape (1,2) or (len(coef), 2).

Returns:

coefs of shifted images.

to_complex(coef)

Return complex valued representation of complex coefficients. This can be useful when comparing, prototyping, or implementing methods from literature.

There is a corresponding method, to_real.

Parameters:

coef – Real Coef from this basis.

Returns:

ComplexCoef representation from this basis.

to_real(complex_coef)

Return real valued representation of complex coefficients. This can be useful when comparing, prototyping, or implementing methods from literature.

There is a corresponding method, to_complex.

Parameters:

complex_coef – Complex Coef from this basis.

Returns:

Real Ceof representation from this basis.

Module contents