GaussianSource

class sofia_redux.scan.source_models.beams.gaussian_source.GaussianSource(peak=1.0, x_mean=0.0, y_mean=0.0, x_fwhm=0.0, y_fwhm=0.0, theta=<Quantity 0. deg>, peak_unit=None, position_unit=None, gaussian_model=None)[source]

Bases: Gaussian2D

Initialize a GaussianSource model.

The GaussianSource is an extension of the Gaussian2D model that allows the model to be fit to (or from) a source map.

Parameters:
peakfloat or units.Quantity, optional

The peak amplitude of the Gaussian.

x_meanfloat or units.Quantity, optional

The position of the peak along the x-axis.

y_meanfloat or units.Quantity, optional

The position of the peak along the y-axis.

x_fwhmfloat or units.Quantity, optional

The Full-Width-Half-Max beam width in the x-direction.

y_fwhmfloat or units.Quantity, optional

The Full-Width-Half-Max beam width in the y-direction.

thetafloat or units.Quantity, optional

The rotation of the beam pertaining to x_fwhm and y_fwhm in relation to the actual (x, y) coordinate axis. If a float value is supplied, it is assumed to be in degrees.

peak_unitunits.Unit or units.Quantity or str, optional

The physical units for the peak amplitude. The default is dimensionless.

position_unitunits.Unit or units.Quantity or str, optional

The physical units of all position based parameters (x_mean, y_mean, x_fwhm, y_fwhm)

gaussian_modelGaussian2D, optional

If supplied, extracts all the above parameters from the supplied model.

Attributes Summary

fwhm_rms

Return the FWHM RMS.

fwhm_significance

Return the FWHM significance.

peak_rms

Return the peak rms.

peak_significance

Return the peak significance.

position

Return the (x, y) position coordinate.

referenced_attributes

Return the names of attributes to be referenced rather than copied.

Methods Summary

convolve_with(psf)

Convolve with a given PSF.

copy()

Return a copy of the Gaussian source.

correct(map2d)

Apply peak value correction from a map.

deconvolve_with(psf)

Deconvolve with a given psf.

edit_header(header[, fits_id, beam_name, ...])

Edit a FITS header with the beam parameters.

find_local_centroid(image)

Find the local peak in a map using the 'centroid' method.

find_local_peak(image[, sign])

Find the local peak in a map using the 'position' method.

find_peak(image[, grid, sign])

Fit the peak coordinates from a given map.

find_source_extent(image[, max_iterations, ...])

Find the extent of the source and shape.

fit_map(map2d[, max_iterations, ...])

Fit the Gaussian source to a given map.

fit_map_least_squares(map2d[, degree, ...])

Fit the Gaussian to a given map using LSQ method (adaptTo).

gaussian_2d_fit(coordinates, amplitude, x0, ...)

A simple 2-dimensional Gaussian function.

get_asymmetry_2d(image, angle, radial_range)

Return the asymmetry of the GaussianSource.

get_center_offset([offset])

Find the offset of the source position from the grid reference.

get_correction_factor(map2d)

Get the correction factor for a given map.

get_data(map2d[, size_unit])

Return a dictionary of properties for to the source model on a map.

get_gaussian_2d()

Return a representation of the beam as a Gaussian2D object.

get_grid_coordinates(index)

Return the grid coordinates for the given map indices.

get_integral(psf_area)

Return the integral over a given area.

get_lsq_fit_parameters(image)

Return the LSQ fit parameters for curve_fit.

get_representation(grid)

Return a representation of the Gaussian source on a new grid.

pointing_info(map2d)

Return a list of strings with pointing information.

scale_fwhm(factor)

Scale the fwhm by a given amount

scale_peak(factor)

Scale the peak by a given amount.

set_center_index(center_index)

Set the center pixel index on the image for the source.

set_exact()

Set the peak to be 'exact'.

set_peak_from(image[, degree, reduce_degrees])

Set the peak value from a given image.

set_peak_position(peak_position)

Set the peak position.

set_positioning_method(method)

Set the peak positioning method.

set_unit(unit)

Set the unit of the Gaussian peak.

uncorrect(map2d)

Uncorrect the peak value for a given map.

Attributes Documentation

fwhm_rms

Return the FWHM RMS.

Returns:
float or astropy.units.Quantity
fwhm_significance

Return the FWHM significance.

Returns:
float
peak_rms

Return the peak rms.

Returns:
float
peak_significance

Return the peak significance.

Returns:
float
position

Return the (x, y) position coordinate.

Returns:
astropy.units.Quantity
referenced_attributes

Return the names of attributes to be referenced rather than copied.

Returns:
set of str

Methods Documentation

convolve_with(psf)[source]

Convolve with a given PSF.

Convolve the current beam with a given PSF. The shape of the beam will be accounted for correctly in both (x, y) dimensions unlike the deconvolve_with method.

Parameters:
psfGaussian2D
Returns:
None
copy()[source]

Return a copy of the Gaussian source.

Returns:
GaussianSource
correct(map2d)[source]

Apply peak value correction from a map.

Parameters:
map2dMap2d
Returns:
None
deconvolve_with(psf)[source]

Deconvolve with a given psf.

Note that this process will result in a circular beam. I.e. one in which the FWHM in the x and y-directions are equivalent.

Parameters:
psfGaussian2D
Returns:
None
edit_header(header, fits_id='', beam_name=None, size_unit=None)[source]

Edit a FITS header with the beam parameters.

Parameters:
headerastropy.io.fits.header.Header

The FITS header to edit.

fits_idstr, optional

Not used.

beam_namestr, optional

The name of the beam. Not implemented for the GaussianSource.

size_unitastropy.units.Unit or Quantity or str, optional

If set, convert the major/minor beam values to this unit before setting in the header.

Returns:
None
static find_local_centroid(image)[source]

Find the local peak in a map using the ‘centroid’ method.

Parameters:
imageFlaggedArray
Returns:
peakCoordinate2D

The (x, y) peak position.

static find_local_peak(image, sign=0)[source]

Find the local peak in a map using the ‘position’ method.

Parameters:
imageFlaggedArray

The image to fit.

signint or float, optional

If positive, fit to positive sources in the map. If negative, fit to negative sources in the map. If zero, fit to source magnitudes in the map. Note, this only takes affect if position is set as the positioning_method.

Returns:
peakCoordinate2D

The (x, y) pixel peak position.

find_peak(image, grid=None, sign=0)[source]

Fit the peak coordinates from a given map. (moveTo)

Parameters:
imageFlaggedArray (float)

The image to fit.

gridSphericalGrid

The grid to fit onto.

signint or float, optional

If positive, fit to positive sources in the map. If negative, fit to negative sources in the map. If zero, fit to source magnitudes in the map. Note, this only takes affect if position is set as the positioning_method.

Returns:
peakCoordinate2D

The (x, y) peak position.

find_source_extent(image, max_iterations=40, radius_increment=0.1, tolerance=0.05)[source]

Find the extent of the source and shape.

Calculates the integral of the source, a circular source radius, and a source mask. The algorithm for determining the source iterates using an expanding circular aperture centered on the currently defined center_index attribute. The iteration will be halted, and parameters updated once either the max_iteration number of iterations has been reached, or:

data_sum <= last_data_sum * (1 + tolerance)

The data_sum is the sum of data values inside the circular aperture, and the last_data_sum is that same sum from the previous iteration. Each iteration the radius of the aperture is increased by:

min(pixel_size, radius_increment * current_radius)

The starting radius is equal to the pixel size in the data.

Parameters:
imageFlaggedArray
max_iterationsint, optional

The maximum number of iterations, each of which increases the search radius by radius_increment.

radius_incrementfloat, optional

The factor by which to increase the search radius between iterations.

tolerancefloat, optional

Halt iterations if the change in data sum is less than 1 + tolerance between iterations.

Returns:
None
fit_map(map2d, max_iterations=40, radius_increment=1.1, tolerance=0.05, degree=3, reduce_degrees=False)[source]

Fit the Gaussian source to a given map.

An Observation2D object or Map2D overlap object may be supplied and fit with a Gaussian source. If an Observation2D is supplied, the Gaussian is fit to the significance image, and thus the FWHM rms may also be determined. Otherwise, all other parameters aside from the FWHM weight will be fit for.

Parameters:
map2dMap2D or Observation2D
max_iterationsint, optional

The maximum number of iterations by which to increase the radius.

radius_incrementfloat, optional

The factor by which to increase the radius for each iteration. The initial radius size is set to min(grid.pixel_size).

tolerancefloat, optional

Stop the iterations if (1-tolerance) * last_sum <= sum <= (1+tolerance).

degreeint, optional

The spline degree used to determine the peak value if necessary.

reduce_degreesbool, optional

If True, allow the spline fit to reduce the number of degrees in cases where there are not enough points available to perform the spline fit of degree. If False, a ValueError will be raised if such a fit fails.

Returns:
data_sumfloat

The sum of the source withing the source radius.

fit_map_least_squares(map2d, degree=3, reduce_degrees=False)[source]

Fit the Gaussian to a given map using LSQ method (adaptTo).

Parameters:
map2dMap2D or Observation2D
degreeint, optional

The spline degree used to fit the map peak value.

reduce_degreesbool, optional

If True, allow the spline fit to reduce the number of degrees in cases where there are not enough points available to perform the spline fit of degree. If False, a ValueError will be raised if such a fit fails.

Returns:
data_sumfloat

The sum of the source withing the source radius.

static gaussian_2d_fit(coordinates, amplitude, x0, y0, sigma)[source]

A simple 2-dimensional Gaussian function.

The return value is:

z = A.exp((-(x - x0)^2) / (2 sigma^2)) +

(-(y - y0)^2) / (2 sigma^2)))

Parameters:
coordinates2-tuple (numpy.ndarray)

The (x, y) coordinates to evaluate.

amplitudefloat

The scaling factor.

x0float

The center of the Gaussian in x.

y0float

The center of the Gaussian in y.

sigmafloat

The Gaussian standard deviation.

Returns:
znumpy.ndarray (float)

The derived function value. Will be the same shape as x or y in the coordinates.

get_asymmetry_2d(image, angle, radial_range)[source]

Return the asymmetry of the GaussianSource.

Parameters:
imageImage2D

The image from which to determine the asymmetry.

angleunits.Quantity

The rotation of the image with respect to the GaussianSource.

radial_rangeRange

The range over which the asymmetry should be calculated.

Returns:
asymmetryAsymmetry2D
get_center_offset(offset=None)[source]

Find the offset of the source position from the grid reference.

Parameters:
offsetCoordinate2D, optional

An optional offset into which to place the results.

Returns:
offsetCoordinate2D
get_correction_factor(map2d)[source]

Get the correction factor for a given map.

Parameters:
map2dMap2D
Returns:
correction_factorfloat

Notes

The original CRUSH would always return a value of 1, since the logic based on the map2d.filter_fwhm returns 1 if it is not NaN, and sets the filtering factor to 0 if it is NaN. This has been remedied here, but does not appear to actually be used during the reduction.

get_data(map2d, size_unit=None)[source]

Return a dictionary of properties for to the source model on a map.

The key values returned are:

  • peak: The fitted peak value

  • dpeak: The fitted peak value RMS

  • peakS2N: The peak signal-to-noise ratio

  • int: The integral of the peak on the map

  • dint: The integral rms of the peak on the map

  • intS2N: The significance of the peak on the map

  • FWHM: The full-width-half maximum of the peak

  • dFWHM: The full-width-half-maximum RMS of the peak

Parameters:
map2dMap2D

The map for which to calculate an integral.

size_unitunits.Unit or str, optional

If set, converts FWHM and dFWHM to size_unit.

Returns:
datadict
get_gaussian_2d()[source]

Return a representation of the beam as a Gaussian2D object.

Returns:
Gaussian2D
get_grid_coordinates(index)[source]

Return the grid coordinates for the given map indices.

Parameters:
indexCoordinate2D

The dimensionless map pixel coordinates.

Returns:
coordinatesCoordinate2D

The physical map coordinates for the given index.

get_integral(psf_area)[source]

Return the integral over a given area.

Parameters:
psf_areaastropy.units.Quantity or float

The area over which to evaluate the integral.

Returns:
integral, weight(float, float) or (Quantity, Quantity)
get_lsq_fit_parameters(image)[source]

Return the LSQ fit parameters for curve_fit.

Determines the entry parameters for a non-linear least squares fit to the Gaussian source. These are the initial guess parameters for the peak amplitude, (x, y) mean, and standard deviation. The bounds for these parameters are also returned. For the GaussianSource, the amplitude and centroid location are unbounded, and a minimum bound on the standard deviation is set, based on the currently set FWHM.

Parameters:
imageObservation2D or SignificanceMap
Returns:
function, initial_values, bounds
get_representation(grid)[source]

Return a representation of the Gaussian source on a new grid.

Parameters:
gridGrid2D

The new grid to apply to a copy of this GaussianSource.

Returns:
GaussianSource

A copy of the Gaussian source on the new grid.

pointing_info(map2d)[source]

Return a list of strings with pointing information.

Parameters:
map2dMap2d
Returns:
infolist (str)
scale_fwhm(factor)[source]

Scale the fwhm by a given amount

Parameters:
factorfloat
Returns:
None
scale_peak(factor)[source]

Scale the peak by a given amount.

Parameters:
factorfloat
Returns:
None
set_center_index(center_index)[source]

Set the center pixel index on the image for the source.

Parameters:
center_indexCoordinate2D

The (x, y) center index.

Returns:
None
set_exact()[source]

Set the peak to be ‘exact’.

Returns:
None
set_peak_from(image, degree=3, reduce_degrees=False)[source]

Set the peak value from a given image.

The peak value is determined from spline interpolation on the image at self.center_index. This will also set the weight of the peak as exact (infinity) or will interpolate from a weight map if available.

Parameters:
imageFlaggedArray or Map2D
degreeint, optional

The spline degree to fit.

reduce_degreesbool, optional

If True, allow the spline fit to reduce the number of degrees in cases where there are not enough points available to perform the spline fit of degree. If False, a ValueError will be raised if such a fit fails.

Returns:
None
set_peak_position(peak_position)[source]

Set the peak position.

Parameters:
peak_positionCoordinate2D or Quantity or numpy.ndarray

The (x, y) peak position.

Returns:
None
set_positioning_method(method)[source]

Set the peak positioning method.

Parameters:
methodstr

May be one of {‘position’, ‘peak’, ‘centroid’}.

Returns:
None
set_unit(unit)[source]

Set the unit of the Gaussian peak.

Parameters:
unitunits.Quantity or units.Unit or float or str or None
Returns:
None
uncorrect(map2d)[source]

Uncorrect the peak value for a given map.

Parameters:
map2dMap2d
Returns:
None