sourcefinder.testutil.convolve#
Functions#
|
Analytically convolve two 2D anisotropic Gaussian profiles. |
|
Evaluate a 2D anisotropic Gaussian defined by its covariance matrix. |
|
From covariance matrix Sigma (2x2) return (sigma_maj, sigma_min, |
Module Contents#
- sourcefinder.testutil.convolve.convolve_gaussians(I1, mu1, Sigma1, I2, mu2, Sigma2)[source]#
Analytically convolve two 2D anisotropic Gaussian profiles.
The convolution of two Gaussians:
G₁(r) = I₁ exp[-½ (r − μ₁)ᵀ Σ₁⁻¹ (r − μ₁)] and G₂(r) = I₂ exp[-½ (r − μ₂)ᵀ Σ₂⁻¹ (r − μ₂)]
yields another Gaussian with parameters:
Σ_H = Σ₁ + Σ₂ μ_H = μ₁ + μ₂ I_H = I₁ I₂ √(|Σ₁| |Σ₂| / |Σ_H|)
- Parameters:
- I1, I2float
Peak brightnesses of the two input Gaussians.
- mu1, mu2array_like, shape (2,)
Center coordinates (x₀, y₀) of each Gaussian.
- Sigma1, Sigma2array_like, shape (2, 2)
Covariance matrices describing the shape and orientation of each Gaussian. Must be symmetric and positive definite.
- Returns:
- I_Hfloat
Peak brightness of the convolved Gaussian.
- mu_Hndarray, shape (2,)
Center coordinates of the convolved Gaussian.
- Sigma_Hndarray, shape (2, 2)
Covariance matrix of the convolved Gaussian.
Notes
The convolution of two Gaussians is itself a Gaussian, which makes this operation fully analytic.
The combined covariance matrix is simply the sum of the individual covariances.
The combined centroid is the vector sum of the input centroids; this assumes both functions are defined in the same coordinate frame (no shift-invariance correction applied).
The normalization factor ensures conservation of integrated flux under convolution.
Examples
>>> import numpy as np >>> I1, mu1, Sigma1 = 1.0, np.array([0., 0.]), np.diag([2.0, 1.0]) >>> I2, mu2, Sigma2 = 1.0, np.array([0.3, -0.2]), np.diag([1.5, 0.5]) >>> I_H, mu_H, Sigma_H = convolve_gaussians(I1, mu1, Sigma1, I2, mu2, Sigma2) >>> I_H, mu_H, Sigma_H (0.707..., array([ 0.3, -0.2]), array([[3.5, 0.], [0., 1.5]]))
- sourcefinder.testutil.convolve.gaussian_from_Sigma_matrix(x, y, I0, mu, Sigma)[source]#
Evaluate a 2D anisotropic Gaussian defined by its covariance matrix.
The Gaussian is given by:
G(x, y) = I₀ * exp[-½ (r - μ)ᵀ Σ⁻¹ (r - μ)],
where Σ is the 2×2 covariance matrix describing the ellipse (its orientation and semi-axes) and μ = (x₀, y₀) is the centroid position.
- Parameters:
- x, yarray_like
Coordinate arrays of identical shape giving the pixel positions at which to evaluate the Gaussian. Typically created with
np.meshgrid.- I0float
Peak intensity (the value of the Gaussian at its center).
- muarray_like, shape (2,)
Center coordinates (x₀, y₀) of the Gaussian, in the same units as x and y.
- Sigmaarray_like, shape (2, 2)
Covariance matrix defining the shape and orientation of the ellipse. Must be symmetric and positive definite.
- Returns:
- Gndarray
The evaluated Gaussian image, with the same shape as x.
Notes
The determinant of Σ controls the spatial extent of the Gaussian.
The integral over all space equals
∫∫ G(x, y) dx dy = I₀ * 2π * sqrt(det(Σ)).
This may be useful for flux normalization.
If Σ is diagonal, the Gaussian is axis-aligned with standard deviations σₓ = sqrt(Σ₀₀) and σᵧ = sqrt(Σ₁₁).
Examples
>>> import numpy as np >>> x, y = np.meshgrid(np.linspace(-5, 5, 101), np.linspace(-5, 5, 101)) >>> mu = np.array([0.0, 0.0]) >>> Sigma = np.array([[4.0, 1.0], ... [1.0, 2.0]]) >>> G = gaussian_from_Sigma_matrix(x, y, I0=1.0, mu=mu, Sigma=Sigma) >>> G.shape (101, 101)
- sourcefinder.testutil.convolve.params_from_sigma(Sigma)[source]#
From covariance matrix Sigma (2x2) return (sigma_maj, sigma_min, theta) with sigma_maj>=sigma_min>=0.
- Parameters:
- Sigmaarray_like, shape (2, 2)
Covariance matrix defining the shape and orientation of the ellipse. Must be symmetric. Covariance matrix elements: if sigma_maj is major axis stddev, sigma_min is minor axis stddev, and theta the position angle (CCW from +Y), then:
S = [[sxx, sxy], [sxy, syy]] = R @ [[sigma_maj^2, 0], [0, sigma_min^2]] @ R.T where R = [[-sin(theta), -cos(theta)], [ cos(theta), -sin(theta)]] i.e. sxx = sigma_maj^2 sin^2(theta) + sigma_min^2 cos^2(theta) syy = sigma_maj^2 cos^2(theta) + sigma_min^2 sin^2(theta) sxy = -(sigma_maj^2 - sigma_min^2) sin(theta) cos(theta)
- Returns:
- sigma_majfloat
Standard deviation along major axis.
- sigma_minfloat
Standard deviation along minor axis.
- thetafloat
Position angle of the major axis in radians, in [0, pi). Measured from the positive y-axis toward the negative-axis.
Notes
The covariance matrix Sigma is related to the ellipse parameters (sigma_maj, sigma_min, theta) by its eigen-decomposition.
The stddevs sigma_maj and sigma_min are the square roots of the eigenvalues of Sigma.
The position angle theta is derived from the eigenvector corresponding to the largest eigenvalue.