AegeanTools modules

angle_tools

Tools for manipulating angles on the surface of a sphere - distance - bearing between two points - translation along a path - paths are either great circles or rhumb lines

also angle <-> string conversion tools for Aegean

AegeanTools.angle_tools.bear(ra1, dec1, ra2, dec2)[source]

Calculate the bearing of point 2 from point 1 along a great circle. The bearing is East of North and is in [0, 360), whereas position angle is also East of North but (-180,180]

Parameters
ra1, dec1, ra2, dec2float

The sky coordinates (degrees) of the two points.

Returns
bearfloat

The bearing of point 2 from point 1 (degrees).

AegeanTools.angle_tools.bear_rhumb(ra1, dec1, ra2, dec2)[source]

Calculate the bearing of point 2 from point 1 along a Rhumb line. The bearing is East of North and is in [0, 360), whereas position angle is also East of North but (-180,180]

Parameters
ra1, dec1, ra2, dec2float

The sky coordinates (degrees) of the two points.

Returns
distfloat

The bearing of point 2 from point 1 along a Rhumb line (degrees).

AegeanTools.angle_tools.dec2dec(dec)[source]

Convert sexegessimal RA string into a float in degrees.

Parameters
decstr

A string separated representing the Dec. Expected format is [+- ]hh:mm[:ss.s] Colons can be replaced with any whit space character.

Returns
decfloat

The Dec in degrees.

AegeanTools.angle_tools.dec2dms(x)[source]

Convert decimal degrees into a sexagessimal string in degrees.

Parameters
xfloat

Angle in degrees

Returns
dmsstr

String of format [+-]DD:MM:SS.SS or XX:XX:XX.XX if x is not finite.

AegeanTools.angle_tools.dec2hms(x)[source]

Convert decimal degrees into a sexagessimal string in hours.

Parameters
xfloat

Angle in degrees

Returns
dmsstring

String of format HH:MM:SS.SS or XX:XX:XX.XX if x is not finite.

AegeanTools.angle_tools.dist_rhumb(ra1, dec1, ra2, dec2)[source]

Calculate the Rhumb line distance between two points [1]. A Rhumb line between two points is one which follows a constant bearing.

Parameters
ra1, dec1, ra2, dec2float

The position of the two points (degrees).

Returns
distfloat

The distance between the two points along a line of constant bearing.

Notes

1

Rhumb line

AegeanTools.angle_tools.gcd(ra1, dec1, ra2, dec2)[source]

Calculate the great circle distance between to points using the haversine formula [1].

Parameters
ra1, dec1, ra2, dec2float

The coordinates of the two points of interest. Units are in degrees.

Returns
distfloat

The distance between the two points in degrees.

Notes

This duplicates the functionality of astropy but is faster as there is no creation of SkyCoords objects.

1

Haversine formula

AegeanTools.angle_tools.ra2dec(ra)[source]

Convert sexegessimal RA string into a float in degrees.

Parameters
rastr

A string separated representing the RA. Expected format is hh:mm[:ss.s] Colons can be replaced with any whit space character.

Returns
rafloat

The RA in degrees.

AegeanTools.angle_tools.translate(ra, dec, r, theta)[source]

Translate a given point a distance r in the (initial) direction theta, along a great circle.

Parameters
ra, decfloat

The initial point of interest (degrees).

r, thetafloat

The distance and initial direction to translate (degrees).

Returns
ra, dec(float, float)

The translated position (degrees).

AegeanTools.angle_tools.translate_rhumb(ra, dec, r, theta)[source]

Translate a given point a distance r in the (initial) direction theta, along a Rhumb line.

Parameters
ra, decfloat

The initial point of interest (degrees).

r, thetafloat

The distance and initial direction to translate (degrees).

Returns
ra, decfloat

The translated position (degrees).

BANE

This module contains all of the BANE specific code The function filter_image should be imported from elsewhere and run as is.

AegeanTools.BANE.barrier(events, sid, kind='neighbour')[source]

act as a multiprocessing barrier

AegeanTools.BANE.filter_image(im_name, out_base, step_size=None, box_size=None, twopass=False, cores=None, mask=True, compressed=False, nslice=None)[source]

Create a background and noise image from an input image. Resulting images are written to outbase_bkg.fits and outbase_rms.fits

Parameters
im_namestr

Image to filter.

out_basestr or None

The output filename base. Will be modified to make _bkg and _rms files. If None, then no files are written.

step_size(int,int)

Tuple of the x,y step size in pixels

box_size(int,int)

The size of the box in pixels

twopassbool

Perform a second pass calculation to ensure that the noise is not contaminated by the background. Default = False

coresint

Number of CPU corse to use. Default = all available

nsliceint

The image will be divided into this many horizontal stripes for processing. Default = None = equal to cores

maskbool

Mask the output array to contain np.nna wherever the input array is nan or not finite. Default = true

compressedbool

Return a compressed version of the background/noise images. Default = False

Returns
bkg, rmsnumpy.ndarray

The computed background and rms maps (not compressed)

AegeanTools.BANE.filter_mc_sharemem(filename, step_size, box_size, cores, shape, nslice=None, domask=True)[source]

Calculate the background and noise images corresponding to the input file. The calculation is done via a box-car approach and uses multiple cores and shared memory.

Parameters
filenamestr

Filename to be filtered.

step_size(int, int)

Step size for the filter.

box_size(int, int)

Box size for the filter.

coresint

Number of cores to use. If None then use all available.

nsliceint

The image will be divided into this many horizontal stripes for processing. Default = None = equal to cores

shape(int, int)

The shape of the image in the given file.

domaskbool

True(Default) = copy data mask to output.

Returns
bkg, rmsnumpy.ndarray

The interpolated background and noise images.

AegeanTools.BANE.get_step_size(header)[source]

Determine the grid spacing for BANE operation.

This is set to being 4x the synthesized beam width. If the beam is not circular then the “width” is sqrt(a*b)

For the standard 4 pix/beam, the step size will be 16 pixels.

Parameters
header
Returns
step_size(int, int)

The grid spacing for BANE operation

AegeanTools.BANE.sigma_filter(filename, region, step_size, box_size, shape, domask, sid)[source]

Calculate the background and rms for a sub region of an image. The results are written to shared memory - irms and ibkg.

Parameters
filenamestring

Fits file to open

regionlist

Region within the fits file that is to be processed. (row_min, row_max).

step_size(int, int)

The filtering step size

box_size(int, int)

The size of the box over which the filter is applied (each step).

shapetuple

The shape of the fits image

domaskbool

If true then copy the data mask to the output.

sidint

The stripe number

Returns
None
AegeanTools.BANE.sigmaclip(arr, lo, hi, reps=3)[source]

Perform sigma clipping on an array, ignoring non finite values.

During each iteration return an array whose elements c obey: mean -std*lo < c < mean + std*hi

where mean/std are the mean std of the input array.

Parameters
arriterable

An iterable array of numeric types.

lofloat

The negative clipping level.

hifloat

The positive clipping level.

repsint

The number of iterations to perform. Default = 3.

Returns
meanfloat

The mean of the array, possibly nan

stdfloat

The std of the array, possibly nan

Notes

Scipy v0.16 now contains a comparable method that will ignore nan/inf values.

AegeanTools.BANE.write_fits(data, header, file_name)[source]

Combine data and a fits header to write a fits file.

Parameters
datanumpy.ndarray

The data to be written.

headerastropy.io.fits.hduheader

The header for the fits file.

file_namestring

The file to write

Returns
None

catalogs

Module for reading at writing catalogs

AegeanTools.catalogs.check_table_formats(files)[source]

Determine whether a list of files are of a recognizable output type.

Parameters
filesstr

A list of file names

Returns
resultbool

True if all the file names are supported

AegeanTools.catalogs.get_table_formats()[source]

Create a list of file extensions that are supported for writing.

Returns
fmtslist

A list of file name extensions that are supported.

AegeanTools.catalogs.load_catalog(filename)[source]

Load a catalogue and extract the source positions (only)

Parameters
filenamestr

Filename to read. Supported types are csv, tab, tex, vo, vot, and xml.

Returns
cataloguelist

A list of [ (ra, dec), …]

AegeanTools.catalogs.load_table(filename)[source]

Load a table from a given file.

Supports csv, tab, tex, vo, vot, xml, fits, and hdf5.

Parameters
filenamestr

File to read

Returns
tableTable

Table of data.

AegeanTools.catalogs.nulls(x)[source]

Convert values of -1 into None.

Parameters
xfloat or int

Value to convert

Returns
val[x, None]
AegeanTools.catalogs.save_catalog(filename, catalog, meta=None, prefix=None)[source]

Save a catalogue of sources using filename as a model. Meta data can be written to some file types (fits, votable).

Each type of source will be in a separate file:

Where filename = base.ext

Parameters
filenamestr

Name of file to write, format is determined by extension.

cataloglist

A list of sources to write. Sources must be of type AegeanTools.models.ComponentSource, AegeanTools.models.SimpleSource, or AegeanTools.models.IslandSource.

prefixstr

Prepend each column name with “prefix_”. Default is to prepend nothing.

metadict

Meta data to be written to the output file. Support for metadata depends on file type.

Returns
None
AegeanTools.catalogs.show_formats()[source]

Print a list of all the file formats that are supported for writing. The file formats are determined by their extensions.

Returns
None
AegeanTools.catalogs.table_to_source_list(table, src_type=<class 'AegeanTools.models.ComponentSource'>)[source]

Convert a table of data into a list of sources.

A single table must have consistent source types given by src_type. src_type should be one of AegeanTools.models.ComponentSource, AegeanTools.models.SimpleSource, or AegeanTools.models.IslandSource.

Parameters
tableTable

Table of sources

src_typeclass

Sources must be of type AegeanTools.models.ComponentSource, AegeanTools.models.SimpleSource, or AegeanTools.models.IslandSource.

Returns
sourceslist

A list of objects of the given type.

AegeanTools.catalogs.update_meta_data(meta=None)[source]

Modify the metadata dictionary. DATE, PROGRAM, and PROGVER are added/modified.

Parameters
metadict

The dictionary to be modified, default = None (empty)

Returns
An updated dictionary.
AegeanTools.catalogs.writeAnn(filename, catalog, fmt)[source]

Write an annotation file that can be read by Kvis (.ann) or DS9 (.reg). Uses ra/dec from catalog. Draws ellipses if bmaj/bmin/pa are in catalog. Draws 30” circles otherwise.

Only AegeanTools.models.ComponentSource will appear in the annotation file unless there are none, in which case AegeanTools.models.SimpleSource (if present) will be written. If any AegeanTools.models.IslandSource objects are present then an island contours file will be written.

Parameters
filenamestr

Output filename base.

cataloglist

List of sources.

fmt[‘ann’, ‘reg’]

Output file type.

Returns
None
AegeanTools.catalogs.writeDB(filename, catalog, meta=None)[source]

Output an sqlite3 database containing one table for each source type

Parameters
filenamestr

Output filename

cataloglist

List of sources of type AegeanTools.models.ComponentSource, AegeanTools.models.SimpleSource, or AegeanTools.models.IslandSource.

metadict

Meta data to be written to table meta

Returns
None
AegeanTools.catalogs.writeFITSTable(filename, table)[source]

Convert a table into a FITSTable and then write to disk.

Parameters
filenamestr

Filename to write.

tableTable

Table to write.

Returns
None

Notes

Due to a bug in numpy, int32 and float32 are converted to int64 and float64 before writing.

AegeanTools.catalogs.writeIslandBoxes(filename, catalog, fmt)[source]

Write an output file in ds9 .reg, or kvis .ann format that contains bounding boxes for all the islands.

Parameters
filenamestr

Filename to write.

cataloglist

List of sources. Only those of type AegeanTools.models.IslandSource will have contours drawn.

fmtstr

Output format type. Currently only ‘reg’ and ‘ann’ are supported. Default = ‘reg’.

Returns
None
AegeanTools.catalogs.writeIslandContours(filename, catalog, fmt='reg')[source]

Write an output file in ds9 .reg format that outlines the boundaries of each island.

Parameters
filenamestr

Filename to write.

cataloglist

List of sources. Only those of type AegeanTools.models.IslandSource will have contours drawn.

fmtstr

Output format type. Currently only ‘reg’ is supported (default)

Returns
None
AegeanTools.catalogs.write_catalog(filename, catalog, fmt=None, meta=None, prefix=None)[source]

Write a catalog (list of sources) to a file with format determined by extension.

Sources must be of type AegeanTools.models.ComponentSource, AegeanTools.models.SimpleSource, or AegeanTools.models.IslandSource.

Parameters
filenamestr

Base name for file to write. _simp, _comp, or _isle will be added to differentiate the different types of sources that are being written.

cataloglist

A list of source objects. Sources must be of type AegeanTools.models.ComponentSource, AegeanTools.models.SimpleSource, or AegeanTools.models.IslandSource.

fmtstr

The file format extension.

prefixstr

Prepend each column name with “prefix_”. Default is to prepend nothing.

metadict

A dictionary to be used as metadata for some file types (fits, VOTable).

Returns
None
AegeanTools.catalogs.write_table(table, filename)[source]

Write a table to a file.

Parameters
tableTable

Table to be written

filenamestr

Destination for saving table.

Returns
None

cluster

Cluster and crossmatch tools and analysis functions.

Includes: - DBSCAN clustering

AegeanTools.cluster.check_attributes_for_regroup(catalog)[source]

Check that the catalog has all the attributes reqired for the regrouping task.

Parameters
cataloglist

List of python objects, ideally derived from AegeanTools.models.SimpleSource

Returns
resultbool

True if the first entry in the catalog has the required attributes

AegeanTools.cluster.norm_dist(src1, src2)[source]

Calculate the normalised distance between two sources. Sources are elliptical Gaussians.

The normalised distance is calculated as the GCD distance between the centers, divided by quadrature sum of the radius of each ellipse along a line joining the two ellipses.

For ellipses that touch at a single point, the normalized distance will be 1/sqrt(2).

Parameters
src1, src2object

The two positions to compare. Objects must have the following parameters: (ra, dec, a, b, pa).

Returns
dist: float

The normalised distance.

AegeanTools.cluster.pairwise_ellpitical_binary(sources, eps, far=None)[source]

Do a pairwise comparison of all sources and determine if they have a normalized distance within eps.

Form this into a matrix of shape NxN.

Parameters
sourceslist

A list of sources (objects with parameters: ra,dec,a,b,pa)

epsfloat

Normalised distance constraint.

farfloat

If sources have a dec that differs by more than this amount then they are considered to be not matched. This is a short-cut around performing GCD calculations.

Returns
probnumpy.ndarray

A 2d array of True/False.

AegeanTools.cluster.regroup(catalog, eps, far=None, dist=<function norm_dist>)[source]

Regroup the islands of a catalog according to their normalised distance. Return a list of island groups. Sources have their (island,source) parameters relabeled.

Parameters
catalogstr or object

Either a filename to read into a source list, or a list of objects with the following properties[units]: ra[deg], dec[deg], a[arcsec], b[arcsec],pa[deg], peak_flux[any]

epsfloat

maximum normalised distance within which sources are considered to be grouped

farfloat

(degrees) sources that are further than this distance appart will not be grouped, and will not be tested. Default = None.

distfunc

a function that calculates the distance between two sources must accept two SimpleSource objects. Default = AegeanTools.cluster.norm_dist()

Returns
islandslist

A list of islands. Each island is a list of sources.

AegeanTools.cluster.regroup_dbscan(srccat, eps=4)[source]

Regroup the islands of a catalog according using DBSCAN.

Return a list of island groups.

Parameters
srccat[ object]

A list of objects with parameters ra,dec (both in decimal degrees)

epsfloat

maximum normalized distance within which sources are considered to be grouped

Returns
islandslist of lists

Each island contains integer indices for members from srccat (in descending dec order).

AegeanTools.cluster.regroup_vectorized(srccat, eps, far=None, dist=<function norm_dist>)[source]

Regroup the islands of a catalog according to their normalised distance.

Assumes srccat is recarray-like for efficiency. Return a list of island groups.

Parameters
srccatnp.rec.arry or pd.DataFrame

Should have the following fields[units]: ra[deg],dec[deg], a[arcsec],b[arcsec],pa[deg], peak_flux[any]

epsfloat

maximum normalised distance within which sources are considered to be grouped

farfloat

(degrees) sources that are further than this distance apart will not be grouped, and will not be tested. Default = 0.5.

distfunc

a function that calculates the distance between a source and each element of an array of sources. Default = AegeanTools.cluster.norm_dist()

Returns
islandslist of lists

Each island contians integer indices for members from srccat (in descending dec order).

AegeanTools.cluster.resize(catalog, ratio=None, psfhelper=None)[source]

Resize all the sources in a given catalogue. Either use a ratio to blindly scale all sources by the same amount, or use a psf map to deconvolve the sources and then convolve them with the new psf

Sources that cannot be rescaled are not returned

Parameters
cataloglist

List of objects

ratiofloat, default=None

Ratio for scaling the sources

psfhelperAegeanTools.wcs_helpers.WCSHelper, default=None

A wcs helper object that contains psf information for the target image/projection

Returns
cataloglist

Modified list of objects

AegeanTools.cluster.sky_dist(src1, src2)[source]

Great circle distance between two sources. A check is made to determine if the two sources are the same object, in this case the distance is zero.

Parameters
src1, src2object

Two sources to check. Objects must have parameters (ra,dec) in degrees.

Returns
distancefloat

The distance between the two sources.

fits_image

Tools for interacting with fits images (HUDLists)

class AegeanTools.fits_image.FitsImage(filename=None, hdu_index=0, beam=None, cube_index=None)[source]

An object that handles the loading and manipulation of a fits file.

get_background_rms()[source]

Calculate the rms of the image. The rms is calculated from the interqurtile range (IQR), to reduce bias from source pixels.

Returns
rmsfloat

The image rms.

Notes

The rms value is cached after first calculation.

get_hdu_header()[source]

Get the image header.

get_pixels()[source]

Get the image data.

Returns
pixelsnumpy.ndarray

2d Array of image pixels.

pix2sky(pixel)[source]

Get the sky coordinates for a given image pixel.

Parameters
pixel(float, float)

Image coordinates.

Returns
ra,decfloat

Sky coordinates (degrees)

set_pixels(pixels)[source]

Set the image data. Will not work if the new image has a different shape than the current image.

Parameters
pixelsnumpy.ndarray

New image data

Returns
None
sky2pix(skypos)[source]

Get the pixel coordinates for a given sky position (degrees).

Parameters
skypos(float,float)

ra,dec position in degrees.

Returns
x,yfloat

Pixel coordinates.

fits_interp

A module to allow fits files to be shrunk in size using decimation, and to be grown in size using interpolation.

AegeanTools.fits_interp.compress(datafile, factor, outfile=None)[source]

Compress a file using decimation.

Parameters
datafilestr or HDUList

Input data to be loaded. (HDUList will be modified if passed).

factorint

Decimation factor.

outfilestr

File to be written. Default = None, which means don’t write a file.

Returns
hdulistHDUList

A decimated HDUList

AegeanTools.fits_interp.expand(datafile, outfile=None)[source]

Expand and interpolate the given data file using the given method. Datafile can be a filename or an HDUList

It is assumed that the file has been compressed and that there are BN_? keywords in the fits header that describe how the compression was done.

Parameters
datafilestr or HDUList

filename or HDUList of file to work on

outfilestr

filename to write to (default = None)

Returns
hdulistHDUList

HDUList of the expanded data.

AegeanTools.fits_interp.load_file_or_hdu(filename)[source]

Load a file from disk and return an HDUList If filename is already an HDUList return that instead

Parameters
filenamestr or HDUList

File or HDU to be loaded

Returns
hdulistHDUList

fitting

Provide fitting routines and helper fucntions to Aegean

AegeanTools.fitting.Bmatrix(C)[source]

Calculate a matrix which is effectively the square root of the correlation matrix C

Parameters
C2d array

A covariance matrix

Returns
B2d array

A matrix B such the B.dot(B’) = inv(C)

AegeanTools.fitting.Cmatrix(x, y, sx, sy, theta)[source]

Construct a correlation matrix corresponding to the data. The matrix assumes a gaussian correlation function.

Parameters
x, yarray-like

locations at which to evaluate the correlation matirx

sx, syfloat

major/minor axes of the gaussian correlation function (sigmas)

thetafloat

position angle of the gaussian correlation function (degrees)

Returns
dataarray-like

The C-matrix.

AegeanTools.fitting.RB_bias(data, pars, ita=None, acf=None)[source]

Calculate the expected bias on each of the parameters in the model pars. Only parameters that are allowed to vary will have a bias. Calculation follows the description of Refrieger & Brown 1998 (cite).

Parameters
data2d-array

data that was fit

parslmfit.Parameters

The model

ita2d-array

The ita matrix (optional).

acf2d-array

The acf for the data.

Returns
biasarray

The bias on each of the parameters

AegeanTools.fitting.bias_correct(params, data, acf=None)[source]

Calculate and apply a bias correction to the given fit parameters

Parameters
paramslmfit.Parameters

The model parameters. These will be modified.

data2d-array

The data which was used in the fitting

acf2d-array

ACF of the data. Default = None.

Returns
None
AegeanTools.fitting.condon_errors(source, theta_n, psf=None)[source]

Calculate the parameter errors for a fitted source using the description of Condon’97 All parameters are assigned errors, assuming that all params were fit. If some params were held fixed then these errors are overestimated.

Parameters
sourceAegeanTools.models.SimpleSource

The source which was fit.

theta_nfloat or None

A measure of the beam sampling. (See Condon’97).

psfAegeanTools.wcs_helpers.Beam

The psf at the location of the source.

Returns
None
AegeanTools.fitting.covar_errors(params, data, errs, B, C=None)[source]

Take a set of parameters that were fit with lmfit, and replace the errors with the 1sigma errors calculated using the covariance matrix.

Parameters
paramslmfit.Parameters

Model

data2d-array

Image data

errs2d-array ?

Image noise.

B2d-array

B matrix.

C2d-array

C matrix. Optional. If supplied then Bmatrix will not be used.

Returns
paramslmfit.Parameters

Modified model.

AegeanTools.fitting.do_lmfit(data, params, B=None, errs=None, dojac=True)[source]

Fit the model to the data data may contain ‘flagged’ or ‘masked’ data with the value of np.NaN

Parameters
data2d-array

Image data

paramslmfit.Parameters

Initial model guess.

B2d-array

B matrix to be used in residual calculations. Default = None.

errs1d-array
dojacbool

If true then an analytic jacobian will be passed to the fitting routine.

Returns
result?

lmfit.minimize result.

paramslmfit.Params

Fitted model.

AegeanTools.fitting.elliptical_gaussian(x, y, amp, xo, yo, sx, sy, theta)[source]

Generate a model 2d Gaussian with the given parameters. Evaluate this model at the given locations x,y.

Parameters
x, ynumeric or array-like

locations at which to evaluate the gaussian

ampfloat

Peak value.

xo, yofloat

Center of the gaussian.

sx, syfloat

major/minor axes in sigmas

thetafloat

position angle (degrees) CCW from x-axis

Returns
datanumeric or array-like

Gaussian function evaluated at the x,y locations.

AegeanTools.fitting.elliptical_gaussian_with_alpha(x, y, v, amp, xo, yo, vo, sx, sy, theta, alpha, beta=None)[source]

Generate a model 2d Gaussian with spectral terms. Evaluate this model at the given locations x,y,dv.

amp is the amplitude at the reference frequency vo

The model is: S(x,v) = amp (v/vo) ^ (alpha + beta x log(v/vo))

When beta is none it is ignored.

Parameters
x, y, vnumeric or array-like

locations at which to evaluate the gaussian

ampfloat

Peak value.

xo, yo, vo: float

Center of the gaussian.

sx, syfloat

major/minor axes in sigmas

thetafloat

position angle (degrees) CCW from x-axis

alpha, beta: float

The spectral terms of the fit.

Returns
datanumeric or array-like

Gaussian function evaluated at the x,y locations.

AegeanTools.fitting.emp_hessian(pars, x, y)[source]

Calculate the hessian matrix empirically. Create a hessian matrix corresponding to the source model ‘pars’ Only parameters that vary will contribute to the hessian. Thus there will be a total of nvar x nvar entries, each of which is a len(x) x len(y) array.

Parameters
parslmfit.Parameters

The model

x, ylist

locations at which to evaluate the Hessian

Returns
hnp.array

Hessian. Shape will be (nvar, nvar, len(x), len(y))

Notes

Uses AegeanTools.fitting.emp_jacobian() to calculate the first order derivatives.

AegeanTools.fitting.emp_jacobian(pars, x, y)[source]

An empirical calculation of the Jacobian Will work for a model that contains multiple Gaussians, and for which some components are not being fit (don’t vary).

Parameters
parslmfit.Model

The model parameters

x, ylist

Locations at which the jacobian is being evaluated

Returns
j2d array

The Jacobian.

AegeanTools.fitting.errors(source, model, wcshelper)[source]

Convert pixel based errors into sky coord errors

Parameters
sourceAegeanTools.models.SimpleSource

The source which was fit.

modellmfit.Parameters

The model which was fit.

wcshelperAegeanTools.wcs_helpers.WCSHelper

WCS information.

Returns
sourceAegeanTools.models.SimpleSource

The modified source obejct.

AegeanTools.fitting.hessian(pars, x, y)[source]

Create a hessian matrix corresponding to the source model ‘pars’ Only parameters that vary will contribute to the hessian. Thus there will be a total of nvar x nvar entries, each of which is a len(x) x len(y) array.

Parameters
parslmfit.Parameters

The model

x, ylist

locations at which to evaluate the Hessian

Returns
hnp.array

Hessian. Shape will be (nvar, nvar, len(x), len(y))

AegeanTools.fitting.jacobian(pars, x, y)[source]

Analytical calculation of the Jacobian for an elliptical gaussian Will work for a model that contains multiple Gaussians, and for which some components are not being fit (don’t vary).

Parameters
parslmfit.Model

The model parameters

x, ylist

Locations at which the jacobian is being evaluated

Returns
j2d array

The Jacobian.

AegeanTools.fitting.lmfit_jacobian(pars, x, y, errs=None, B=None, emp=False)[source]

Wrapper around AegeanTools.fitting.jacobian() and AegeanTools.fitting.emp_jacobian() which gives the output in a format that is required for lmfit.

Parameters
parslmfit.Model

The model parameters

x, ylist

Locations at which the jacobian is being evaluated

errslist

a vector of 1sigma errors (optional). Default = None

B2d-array

a B-matrix (optional) see AegeanTools.fitting.Bmatrix()

empbool

If true the use the empirical Jacobian, otherwise use the analytical one. Default = False.

Returns
j2d-array

A Jacobian.

AegeanTools.fitting.make_ita(noise, acf=None)[source]

Create the matrix ita of the noise where the noise may be a masked array where ita(x,y) is the correlation between pixel pairs that have the same separation as x and y.

Parameters
noise2d-array

The noise image

acf2d-array

The autocorrelation matrix. (None = calculate from data). Default = None.

Returns
ita2d-array

The matrix ita

AegeanTools.fitting.nan_acf(noise)[source]

Calculate the autocorrelation function of the noise where the noise is a 2d array that may contain nans

Parameters
noise2d-array

Noise image.

Returns
acf2d-array

The ACF.

AegeanTools.fitting.new_errors(source, model, wcshelper)[source]

Convert pixel based errors into sky coord errors Uses covariance matrix for ra/dec errors and calculus approach to a/b/pa errors

Parameters
sourceAegeanTools.models.SimpleSource

The source which was fit.

modellmfit.Parameters

The model which was fit.

wcshelperAegeanTools.wcs_helpers.WCSHelper

WCS information.

Returns
sourceAegeanTools.models.SimpleSource

The modified source obejct.

AegeanTools.fitting.ntwodgaussian_lmfit(params)[source]

Convert an lmfit.Parameters object into a function which calculates the model.

Parameters
paramslmfit.Parameters

Model parameters, can have multiple components.

Returns
modelfunc

A function f(x,y) that will compute the model.

flags

Flag constants for use by Aegean.

MIMAS

MIMAS - The Multi-resolution Image Mask for Aegean Software

TODO: Write an in/out reader for MOC formats described by http://arxiv.org/abs/1505.02937

class AegeanTools.MIMAS.Dummy(maxdepth=8)[source]

A state storage class for MIMAS to work with.

Attributes
add_regionlist

List of AegeanTools.MIMAS.Region to be added.

rem_regionlist

List of AegeanTools.MIMAS.Region to be subtracted.

include_circles[[ra, dec, radius],…]

List of circles to be added to the region, units are degrees.

exclude_circles[[ra, dec, radius], …]

List of circles to be subtracted from the region, units are degrees.

include_polygons[[ra,dec, …], …]

List of polygons to be added to the region, units are degrees.

exclude_polygons[[ra,dec, …], …]

List of polygons to be subtracted from the region, units are degrees.

maxdepthint

Depth or resolution of the region for HEALPix. There are 4*2**maxdepth pixels at the deepest layer. Default = 8.

galactic: bool

If true then all ra/dec coordinates will be interpreted as if they were in galactic lat/lon (degrees)

AegeanTools.MIMAS.box2poly(line)[source]

Convert a string that describes a box in ds9 format, into a polygon that is given by the corners of the box

Parameters
linestr

A string containing a DS9 region command for a box.

Returns
poly[ra, dec, …]

The corners of the box in clockwise order from top left.

AegeanTools.MIMAS.circle2circle(line)[source]

Parse a string that describes a circle in ds9 format.

Parameters
linestr

A string containing a DS9 region command for a circle.

Returns
circle[ra, dec, radius]

The center and radius of the circle.

AegeanTools.MIMAS.combine_regions(container)[source]

Return a region that is the combination of those specified in the container. The container is typically a results instance that comes from argparse.

Order of construction is: add regions, subtract regions, add circles, subtract circles, add polygons, subtract polygons.

Parameters
containerAegeanTools.MIMAS.Dummy

The regions to be combined.

Returns
regionAegeanTools.regions.Region

The constructed region.

AegeanTools.MIMAS.galactic2fk5(l, b)[source]

Convert galactic l/b to fk5 ra/dec

Parameters
l, bfloat

Galactic coordinates in radians.

Returns
ra, decfloat

FK5 ecliptic coordinates in radians.

AegeanTools.MIMAS.intersect_regions(flist)[source]

Construct a region which is the intersection of all regions described in the given list of file names.

Parameters
flistlist

A list of region filenames.

Returns
regionAegeanTools.regions.Region

The intersection of all regions, possibly empty.

AegeanTools.MIMAS.mask2mim(maskfile, mimfile, threshold=1.0, maxdepth=8)[source]

Use a fits file as a mask to create a region file.

Pixels in mask file that are equal or above the threshold will be included in the reigon, while those that are below the threshold will not.

Parameters
maskfilestr

Input file in fits format.

mimfilestr

Output filename

thresholdfloat

threshold value for separating include/exclude values

maxdepthint

Maximum depth (resolution) of the healpix pixels

AegeanTools.MIMAS.mask_catalog(regionfile, infile, outfile, negate=False, racol='ra', deccol='dec')[source]

Apply a region file as a mask to a catalog, removing all the rows with ra/dec inside the region If negate=False then remove the rows with ra/dec outside the region.

Parameters
regionfilestr

A file which can be loaded as a AegeanTools.regions.Region. The catalogue will be masked according to this region.

infilestr

Input catalogue.

outfilestr

Output catalogue.

negatebool

If True then pixels outside the region are masked. Default = False.

racol, deccolstr

The name of the columns in table that should be interpreted as ra and dec. Default = ‘ra’, ‘dec’

AegeanTools.MIMAS.mask_file(regionfile, infile, outfile, negate=False)[source]

Created a masked version of file, using a region.

Parameters
regionfilestr

A file which can be loaded as a AegeanTools.regions.Region. The image will be masked according to this region.

infilestr

Input FITS image.

outfilestr

Output FITS image.

negatebool

If True then pixels outside the region are masked. Default = False.

AegeanTools.MIMAS.mask_plane(data, wcs, region, negate=False)[source]

Mask a 2d image (data) such that pixels within ‘region’ are set to nan.

Parameters
data2d-array

Image array.

wcsastropy.wcs.WCS

WCS for the image in question.

regionAegeanTools.regions.Region

A region within which the image pixels will be masked.

negatebool

If True then pixels outside the region are masked. Default = False.

Returns
masked2d-array

The original array, but masked as required.

AegeanTools.MIMAS.mask_table(region, table, negate=False, racol='ra', deccol='dec')[source]

Apply a given mask (region) to the table, removing all the rows with ra/dec inside the region If negate=False then remove the rows with ra/dec outside the region.

Parameters
regionAegeanTools.regions.Region

Region to mask.

tableAstropy.table.Table

Table to be masked.

negatebool

If True then pixels outside the region are masked. Default = False.

racol, deccolstr

The name of the columns in table that should be interpreted as ra and dec. Default = ‘ra’, ‘dec’

Returns
maskedAstropy.table.Table

A view of the given table which has been masked.

AegeanTools.MIMAS.mim2fits(mimfile, fitsfile)[source]

Convert a MIMAS region (.mim) file into a MOC region (.fits) file.

Parameters
mimfilestr

Input file in MIMAS format.

fitsfilestr

Output file.

AegeanTools.MIMAS.mim2reg(mimfile, regfile)[source]

Convert a MIMAS region (.mim) file into a DS9 region (.reg) file.

Parameters
mimfilestr

Input file in MIMAS format.

regfilestr

Output file.

AegeanTools.MIMAS.poly2poly(line)[source]

Parse a string of text containing a DS9 description of a polygon.

This function works but is not very robust due to the constraints of healpy.

Parameters
linestr

A string containing a DS9 region command for a polygon.

Returns
poly[ra, dec, …]

The coordinates of the polygon.

AegeanTools.MIMAS.reg2mim(regfile, mimfile, maxdepth)[source]

Parse a DS9 region file and write a MIMAS region (.mim) file.

Parameters
regfilestr

DS9 region (.reg) file.

mimfilestr

MIMAS region (.mim) file.

maxdepthstr

Depth/resolution of the region file.

AegeanTools.MIMAS.save_as_image(region, filename)[source]

Convert a MIMAS region (.mim) file into a image (eg .png)

Parameters
regionAegeanTools.regions.Region

Region of interest.

filenamestr

Output filename.

AegeanTools.MIMAS.save_region(region, filename)[source]

Save the given region to a file

Parameters
regionAegeanTools.regions.Region

A region.

filenamestr

Output file name.

models

Different types of sources that Aegean is able to fit

class AegeanTools.models.ComponentSource[source]

A Gaussian component, aka a source, that was measured by Aegean.

Attributes
islandint

The island which this component is part of.

sourceint

The source number within the island.

background, local_rmsfloat

Background and local noise level in the image at the location of this source.

ra, err_ra, dec, err-decfloat

Sky location of the source including uncertainties. Decimal degrees.

ra_str, dec_strstr

Sky location in HH:MM:SS.SS +DD:MM:SS.SS format.

galacticbool

If true then ra,dec are interpreted as glat,glon instead. Default = False. This is a class attribute, not an instance attribute.

peak_flux, err_peak_fluxfloat

The peak flux and associated uncertainty.

int_flux, err_int_fluxfloat

Integrated flux and associated uncertainty.

a, err_a, b, err_b, pa, err_pa: float

Shape parameters for this source and associated uncertainties. a/b are in arcsec, pa is in degrees East of North.

residual_mean, residual_stdfloat

The mean and standard deviation of the model-data for this island of pixels.

psf_a, psf_b, psf_pafloat

The shape parameters for the point spread function (degrees).

flagsint

Flags. See AegeanTools.flags.

uuidstr

Unique ID for this source. This is random and not dependent on the source properties.

class AegeanTools.models.DummyLM[source]

A dummy copy of the lmfit results, for use when no fitting was done.

Attributes
residual[np.nan, np.nan]

The residual background and rms.

success: bool

False - the fitting has failed.

class AegeanTools.models.GlobalFittingData[source]

A class to hold the properties associated with an image. [ These were once in the global scope of a monolithic script, hence the name]. (should be) Read-only once created. Used by island fitting subprocesses.

Attributes
imgAegeanTools.fits_image.FitsImage

Image that is being analysed, aka the input image.

dcurve2d-array

Image of +1,0,-1 representing the curvature of the input image.

rmsimg, bkgimg2d-array

The noise and background of the input image.

hdu_headerHDUHeader

FITS header for the input image.

beamAegeanTools.wcs_helpers.Beam

The synthesized beam of the input image.

data_pix2d-array

A link to the data array that is contained within the img.

dtype{np.float32, np.float64}

The data type for the input image. Will be enforced upon writing.

regionAegeanTools.regions.Region

The region that will be used to limit the source finding of Aegean.

wcshelperAegeanTools.wcs_helpers.WCSHelper

A helper object for WCS operations, created from hdu_header.

blankbool

If true, then the input image will be blanked at the location of each of the measured islands.

class AegeanTools.models.IslandFittingData(isle_num=0, i=None, scalars=None, offsets=(0, 0, 1, 1), doislandflux=False)[source]

All the data required to fit a single island. Instances are pickled and passed to the fitting subprocesses

Attributes
isle_numint

island number

i2d-array

a 2D numpy array of pixel values

scalars(innerclip, outerclip, max_summits)

Inner and outer clipping limits (sigma), and the maximum number of components that should be fit.

offsets(xmin, xmax, ymin, ymax)

The offset between the boundaries of the island i, within the larger image.

doislandfluxboolean

If true then also measure properties of the island.

class AegeanTools.models.IslandSource[source]

An island of pixels.

Attributes
island: int

The island identification number

componentsint

The number of components that make up this island.

background, local_rmsfloat

Background and local noise level in the image at the location of this source.

ra, decfloat

Sky location of the brightest pixel in this island. Decimal degrees.

ra_str, dec_strstr

Sky location in HH:MM:SS.SS +DD:MM:SS.SS format.

galacticbool

If true then ra,dec are interpreted as glat,glon instead. Default = False. This is a class attribute, not an instance attribute.

peak_flux, peak_pixelfloat

Value of the brightest pixel for this source.

int_flux, err_int_fluxfloat

Integrated flux and associated uncertainty.

x_width, y_widthint

The extent of the island in pixel space. The width is of the smallest bounding box.

max_angular_sizefloat

The maximum angular size of the island in sky coordinates (degrees).

pafloat

Position angle for the line representing the maximum angular size.

pixelsint

The number of pixels covered by this island.

areafloat

The area of this island in sky coordinates (square degrees).

beam_areafloat

The area of the synthesized beam of the image at the location of the brightest pixel. (square degrees).

etafloat

A factor that accounts for the difference between the integrated flux counted by summing pixels, and the integrated flux that would be produced by integrating an appropriately sized Gaussian.

extentfloat
contourlist

A list of pixel coordinates that mark the pixel boundaries for this island of pixels.

max_angular_size_anchors[x1, y1, x2, y2]

The end points of the vector that describes the maximum angular size of this island.

flagsint

Flags. See AegeanTools.flags.

uuidstr

Unique ID for this source. This is random and not dependent on the source properties.

class AegeanTools.models.PixelIsland(dim=2)[source]

An island of pixels within an image or cube

Attributes
dimint

The number of dimensions of this island. dim >=2, default is 2 (ra/dec).

bounding_box[(min, max), (min, max), …]

A bounding box for this island. len(bounding_box)==dim.

masknp.array(dtype=bool)

A mask that represents the island within the bounding box.

calc_bounding_box(data, offsets)[source]

Compute the bounding box for a data cube of dimension dim. The bounding box will be the smallest nd-cube that bounds the non-zero entries of the cube.

Parameters
datanp.ndarray

Data array with dimension equal to self.dim

offsets[xmin, ymin, …]

The offset between the image zero index and the zero index of data. len(offsets)==dim

set_mask(data)[source]
Parameters
datanp.array
class AegeanTools.models.SimpleSource[source]

The base source class for an elliptical Gaussian.

Attributes
background, local_rmsfloat

Background and local noise level in the image at the location of this source.

ra, decfloat

Sky location of this source. Decimal degrees.

galacticbool

If true then ra,dec are interpreted as glat,glon instead. Default = False. This is a class attribute, not an instance attribute.

peak_flux, err_peak_fluxfloat

The peak flux value and associated uncertainty.

peak_pixelfloat

Value of the brightest pixel for this source.

flagsint

Flags. See AegeanTools.flags.

a, b, pafloat

Shape parameters for this source.

uuidstr

Unique ID for this source. This is random and not dependent on the source properties.

as_list()[source]

Return an ordered list of the source attributes

AegeanTools.models.classify_catalog(catalog)[source]

Look at a list of sources and split them according to their class.

Parameters
catalogiterable

A list or iterable object of {SimpleSource, IslandSource, ComponentSource} objects, possibly mixed. Any other objects will be silently ignored.

Returns
componentslist

List of sources of type ComponentSource

islandslist

List of sources of type IslandSource

simpleslist

List of source of type SimpleSource

AegeanTools.models.island_itergen(catalog)[source]

Iterate over a catalog of sources, and return an island worth of sources at a time. Yields a list of components, one island at a time

Parameters
catalogiterable

A list or iterable of AegeanTools.models.ComponentSource objects.

Yields
grouplist

A list of all sources within an island, one island at a time.

msq2

Provie a class which performs the marching squares algorithm on an image. The desired output is a set of regions / contours.

class AegeanTools.msq2.MarchingSquares(data)[source]

Implementation of a marching squares algorithm. With reference to http://devblog.phillipspiess.com/2010/02/23/better-know-an-algorithm-1-marching-squares/ but written in python

do_march()[source]

March about and trace the outline of our object

Returns
perimeterlist

The pixels on the perimeter of the region [[x1, y1], …]

do_march_all()[source]

Recursive march in the case that we have a fragmented shape.

Returns
perimeters[perimeter1, …]

The perimeters of all the regions in the image.

find_start_point()[source]

Find the first location in our array that is not empty

solid(x, y)[source]

Determine whether the pixel x,y is nonzero

Parameters
x, yint

The pixel of interest.

Returns
solidbool

True if the pixel is not zero.

step(x, y)[source]

Move from the current location to the next

Parameters
x, yint

The current location

walk_perimeter(startx, starty)[source]

Starting at a point on the perimeter of a region, ‘walk’ the perimeter to return to the starting point. Record the path taken.

Parameters
startx, startyint

The starting location. Assumed to be on the perimeter of a region.

Returns
perimeterlist

A list of pixel coordinates [ [x1,y1], …] that constitute the perimeter of the region.

regions

Describe sky areas as a collection of HEALPix pixels

class AegeanTools.regions.Region(maxdepth=11)[source]

A Region object represents a footprint on the sky. This is done in a way similar to a MOC. The region is stored as a list of healpix pixels, allowing for binary set-like operations.

Attributes
maxdepthint

The depth or resolution of the region. At the deepest level there will be 4*2**maxdepth pixels on the sky. Default = 11

pixeldictdict

A dictionary of sets, each set containing the pixels within the region. The sets are indexed by their layer number.

demotedset

A representation of this region at the deepest layer.

add_circles(ra_cen, dec_cen, radius, depth=None)[source]

Add one or more circles to this region

Parameters
ra_cen, dec_cen, radiusfloat or list

The center and radius of the circle or circles to add to this region.

depthint

The depth at which the given circles will be inserted.

add_pixels(pix, depth)[source]

Add one or more HEALPix pixels to this region.

Parameters
pixint or iterable

The pixels to be added

depthint

The depth at which the pixels are added.

add_poly(positions, depth=None)[source]

Add a single polygon to this region.

Parameters
positions[[ra, dec], …]

Positions for the vertices of the polygon. The polygon needs to be convex and non-intersecting.

depthint

The deepth at which the polygon will be inserted.

get_area(degrees=True)[source]

Calculate the total area represented by this region.

Parameters
degreesbool

If True then return the area in square degrees, otherwise use steradians. Default = True.

Returns
areafloat

The area of the region.

get_demoted()[source]

Get a representation of this region at the deepest level.

Returns
demotedset

A set of pixels, at the highest resolution.

intersect(other)[source]

Combine with another Region by performing intersection on their pixlists.

Requires both regions to have the same maxdepth.

Parameters
otherAegeanTools.regions.Region

The region to be combined.

classmethod load(mimfile)[source]

Create a region object from the given file.

Parameters
mimfilestr

File to load.

Returns
regionAegeanTools.regions.Region

A region object

static radec2sky(ra, dec)[source]

Convert [ra], [dec] to [(ra[0], dec[0]),….] and also ra,dec to [(ra,dec)] if ra/dec are not iterable

Parameters
ra, decfloat or iterable

Sky coordinates

Returns
skynumpy.array

array of (ra,dec) coordinates.

save(mimfile)[source]

Save this region to a file

Parameters
mimfilestr

File to write

static sky2ang(sky)[source]

Convert ra,dec coordinates to theta,phi coordinates ra -> phi dec -> theta

Parameters
skynumpy.array

Array of (ra,dec) coordinates. See AegeanTools.regions.Region.radec2sky()

Returns
theta_phinumpy.array

Array of (theta,phi) coordinates.

classmethod sky2vec(sky)[source]

Convert sky positions in to 3d-vectors on the unit sphere.

Parameters
skynumpy.array

Sky coordinates as an array of (ra,dec)

Returns
vecnumpy.array

Unit vectors as an array of (x,y,z)

sky_within(ra, dec, degin=False)[source]

Test whether a sky position is within this region

Parameters
ra, decfloat

Sky position.

deginbool

If True the ra/dec is interpreted as degrees, otherwise as radians. Default = False.

Returns
withinbool

True if the given position is within one of the region’s pixels.

symmetric_difference(other)[source]

Combine with another Region by performing the symmetric difference of their pixlists.

Requires both regions to have the same maxdepth.

Parameters
otherAegeanTools.regions.Region

The region to be combined.

union(other, renorm=True)[source]

Add another Region by performing union on their pixlists.

Parameters
otherAegeanTools.regions.Region

The region to be combined.

renormbool

Perform renormalisation after the operation? Default = True.

classmethod vec2sky(vec, degrees=False)[source]

Convert [x,y,z] vectors into sky coordinates ra,dec

Parameters
vecnumpy.array

Unit vectors as an array of (x,y,z)

degrees
Returns
skynumpy.array

Sky coordinates as an array of (ra,dec)

without(other)[source]

Subtract another Region by performing a difference operation on their pixlists.

Requires both regions to have the same maxdepth.

Parameters
otherAegeanTools.regions.Region

The region to be combined.

write_fits(filename, moctool='')[source]

Write a fits file representing the MOC of this region.

Parameters
filenamestr

File to write

moctoolstr

String to be written to fits header with key “MOCTOOL”. Default = ‘’

write_reg(filename)[source]

Write a ds9 region file that represents this region as a set of diamonds.

Parameters
filenamestr

File to write

source_finder

The Aegean source finding program.

class AegeanTools.source_finder.SourceFinder(**kwargs)[source]

The Aegean source finding algorithm

Attributes
global_dataAegeanTools.models.GlobalFittingData

State holder for properties.

sourceslist

List of sources that have been found/measured.

loglogging.log

Logger to use. Default = None

estimate_lmfit_parinfo(data, rmsimg, curve, beam, innerclip, outerclip=None, offsets=(0, 0), max_summits=None)[source]

Estimates the number of sources in an island and returns initial parameters for the fit as well as limits on those parameters.

Parameters
data2d-array

(sub) image of flux values. Background should be subtracted.

rmsimg2d-array

Image of 1sigma values

curve2d-array

Image of curvature values [-1,0,+1]

beamAegeanTools.fits_image.Beam

The beam information for the image.

innerclip, outerclipfloat

Inerr and outer level for clipping (sigmas).

offsets(int, int)

The (x,y) offset of data within it’s parent image

max_summitsint

If not None, only this many summits/components will be fit. More components may be present in the island, but subsequent components will not have free parameters.

Returns
modellmfit.Parameters

The initial estimate of parameters for the components within this island.

find_sources_in_image(filename, hdu_index=0, outfile=None, rms=None, bkg=None, max_summits=None, innerclip=5, outerclip=4, cores=None, rmsin=None, bkgin=None, beam=None, doislandflux=False, nopositive=False, nonegative=False, mask=None, imgpsf=None, blank=False, docov=True, cube_index=None, progress=False)[source]

Run the Aegean source finder.

Parameters
filenamestr or HDUList

Image filename or HDUList.

hdu_indexint

The index of the FITS HDU (extension).

outfilestr

file for printing catalog (NOT a table, just a text file of my own design)

rmsfloat

Use this rms for the entire image (will also assume that background is 0)

max_summitsint

Fit up to this many components to each island (extras are included but not fit)

innerclip, outerclipfloat

The seed (inner) and flood (outer) clipping level (sigmas).

coresint

Number of CPU cores to use. None means all cores.

rmsin, bkginstr or HDUList

Filename or HDUList for the noise and background images. If either are None, then it will be calculated internally.

beam(major, minor, pa)

Floats representing the synthesised beam (degrees). Replaces whatever is given in the FITS header. If the FITS header has no BMAJ/BMIN then this is required.

doislandfluxbool

If True then each island will also be characterized.

nopositive, nonegativebool

Whether to return positive or negative sources. Default nopositive=False, nonegative=True.

maskstr

The filename of a region file created by MIMAS. Islands outside of this region will be ignored.

imgpsfstr or HDUList

Filename or HDUList for a psf image.

blankbool

Cause the output image to be blanked where islands are found.

docovbool

If True then include covariance matrix in the fitting process. Default=True

cube_indexint

For image cubes, cube_index determines which slice is used.

progressbool

If true then show a progress bar when fitting island groups

Returns
sourceslist

List of sources found.

load_globals(filename, hdu_index=0, bkgin=None, rmsin=None, beam=None, verb=False, rms=None, bkg=None, cores=1, do_curve=False, mask=None, psf=None, blank=False, docov=True, cube_index=None)[source]

Populate the global_data object by loading or calculating the various components

Parameters
filenamestr or HDUList

Main image which source finding is run on

hdu_indexint

HDU index of the image within the fits file, default is 0 (first)

bkgin, rmsinstr or HDUList

background and noise image filename or HDUList

beamAegeanTools.fits_image.Beam

Beam object representing the synthsized beam. Will replace what is in the FITS header.

verbbool

Verbose. Write extra lines to INFO level log.

rms, bkgfloat

A float that represents a constant rms/bkg levels for the image. Default = None, which causes the rms/bkg to be loaded or calculated.

coresint

Number of cores to use if different from what is autodetected.

do_curvebool

If True a curvature map will be created, default=True.

maskstr or AegeanTools.regions.Region

filename or Region object

psfstr or HDUList

Filename or HDUList of a psf image

blankbool

True = blank output image where islands are found. Default = False.

docovbool

True = use covariance matrix in fitting. Default = True.

cube_indexint

For an image cube, which slice to use.

priorized_fit_islands(filename, catalogue, hdu_index=0, outfile=None, bkgin=None, rmsin=None, cores=1, rms=None, bkg=None, beam=None, imgpsf=None, catpsf=None, stage=3, ratio=None, outerclip=3, doregroup=True, regroup_eps=None, docov=True, cube_index=None, progress=False)[source]

Take an input catalog, and image, and optional background/noise images fit the flux and ra/dec for each of the given sources, keeping the morphology fixed

Multiple cores can be specified, and will be used.

Parameters
filenamestr or HDUList

Image filename or HDUList.

cataloguestr or list

Input catalogue file name or list of ComponentSource objects.

hdu_indexint

The index of the FITS HDU (extension).

outfilestr

file for printing catalog (NOT a table, just a text file of my own design)

rmsin, bkginstr or HDUList

Filename or HDUList for the noise and background images. If either are None, then it will be calculated internally.

coresint

Number of CPU cores to use. None means all cores.

rmsfloat

Use this rms for the entire image (will also assume that background is 0)

beam(float, float, float)

(major, minor, pa) representing the synthesised beam (degrees). Replaces whatever is given in the FITS header. If the FITS header has no BMAJ/BMIN then this is required.

imgpsfstr or HDUList

Filename or HDUList for a psf image.

catpsfstr or HDUList

Filename or HDUList for the catalogue psf image.

stageint

Refitting stage

ratiofloat

If not None - ratio of image psf to catalog psf, otherwise interpret from catalogue or image if possible

outerclipfloat

The flood (outer) clipping level (sigmas).

doregroupbool, Default=True

Relabel all the islands/groups to ensure that nearby components are jointly fit.

regroup_eps: float, Default=None

The linking parameter for regouping. Components that are closer than this distance (in arcmin) will be jointly fit. If NONE, then use 4x the average source major axis size (after rescaling if required).

docovbool, Default=True

If True then include covariance matrix in the fitting process.

cube_indexint, Default=None

For image cubes, slice determines which slice is used.

progressbool, Default=True

Show a progress bar when fitting island groups

Returns
sourceslist

List of sources measured.

result_to_components(result, model, island_data, isflags)[source]

Convert fitting results into a set of components

Parameters
resultlmfit.MinimizerResult

The fitting results.

modellmfit.Parameters

The model that was fit.

island_dataAegeanTools.models.IslandFittingData

Data about the island that was fit.

isflagsint

Flags that should be added to this island in addition to those within the model)

Returns
sourceslist

A list of components, and islands if requested.

save_background_files(image_filename, hdu_index=0, bkgin=None, rmsin=None, beam=None, rms=None, bkg=None, cores=1, outbase=None)[source]

Generate and save the background and RMS maps as FITS files. They are saved in the current directly as aegean-background.fits and aegean-rms.fits.

Parameters
image_filenamestr or HDUList

Input image.

hdu_indexint

If fits file has more than one hdu, it can be specified here. Default = 0.

bkgin, rmsinstr or HDUList

Background and noise image filename or HDUList

beamAegeanTools.fits_image.Beam

Beam object representing the synthsized beam. Will replace what is in the FITS header.

rms, bkgfloat

A float that represents a constant rms/bkg level for the image. Default = None, which causes the rms/bkg to be loaded or calculated.

coresint

Number of cores to use if different from what is autodetected.

outbasestr

Basename for output files.

save_image(outname)[source]

Save the image data. This is probably only useful if the image data has been blanked.

Parameters
outnamestr

Name for the output file.

AegeanTools.source_finder.characterise_islands(islands, im, bkg, rms, wcshelper, err_type='best', max_summits=None, do_islandfit=False)[source]

Do the source characterisation based on the initial estimate of the island properties.

Parameters
islands[lmfit.Parameters, … ]

The initial estimate of parameters for the components within each island.

im, bkg, rmsnp.ndarray

The image, background, and noise maps

wcshelperAegeanTools.wcs_helpers.WCSHelper

A wcs helper object

err_typestr or None

The method for calculating uncertainties on parameters: ‘best’ - Uncertainties measured based on covariance matrix of the

fit and of the data See Hancock et al. 2018 for a description of this process.

‘condon’ - Uncertainties are calculated based on Condon’98 (?year) ‘raw’ - uncertainties directly from the covariance matrix only ‘none’ or None - No uncertainties, all will be set to -1.

max_summitsint

The maximum number of summits that will be fit. The final model may contain additional components but only the first few will be fit.

do_islandfitbool

If True, then also characterise islands as well as components. Default=False.

Returns
sources[AegeanTools.models.SimpleSource, … ]

A list of characterised sources of type AegeanTools.models.impleSource, AegeanTools.models.ComponentSource, or AegeanTools.models.IslandSource.

AegeanTools.source_finder.check_cores(cores)[source]

Determine how many cores we are able to use. Return 1 if we are not able to make a queue via pprocess.

Parameters
coresint

The number of cores that are requested.

Returns
coresint

The number of cores available.

AegeanTools.source_finder.estimate_parinfo_image(islands, im, rms, wcshelper, max_summits=None, log=<Logger dummy (WARNING)>)[source]

Estimate the initial parameters for fitting for each of the islands of pixels. The source sizes will be initialised as the psf of the image, which is either determined by the WCS of the image file or the psf map if one is supplied.

Parameters
islands[AegeanTools.models.IslandFittingData, … ]

A list of islands which will be converted into groups of sources

im, rmsnumpy.ndarray

The image and noise maps

wcshelperAegeanTools.wcs_helpers.WCSHelper

A wcshelper object valid for the image map

max_summitsint or None

The maximum number of summits that will be fit. Any in addition to this will be estimated but their parameters will have vary=False.

loglogging.Logger or None

For handling logs (or not)

max_summitsint

The maximum number of summits that will be fit. Any in addition to this will be estimated but their parameters will have vary=False.

loglogging.Logger or None

For handling logs (or not)

Returns
sources[lmfit.Parameters, … ]

The initial estimate of parameters for the components within each island.

AegeanTools.source_finder.find_islands(im, bkg, rms, seed_clip=5.0, flood_clip=4.0, log=<Logger dummy (WARNING)>)[source]

This function designed to be run as a stand alone process

Parameters
im, bkg, rmsnumpy.ndarray

Image, background, and rms maps

seed_clip, flood_clipfloat

The seed clip which is used to create islands, and flood clip which is used to grow islands. The units are in SNR.

loglogging.Logger or None

For handling logs (or not)

Returns
islands[AegeanTools.models.PixelIsland, …]

a list of islands

AegeanTools.source_finder.fit_islands_parinfo(models, im, rms, wcshelper)[source]

Turn a list of sources into a set of islands and parameter estimates which can then be characterised.

Parameters
models[lmfit.Parinfo, … ]

A list of sources in the catalogue.

imnp.ndarray

The image map

wcshelperAegeanTools.wcs_helpers.WCSHelper

A wcs object valid for the image map

Returns
islands[AegeanTools.models.SimpleSource, …]

a list of islands

AegeanTools.source_finder.fix_shape(source)[source]

Ensure that a>=b for a given source object. If a<b then swap a/b and increment pa by 90. err_a/err_b are also swapped as needed.

Parameters
sourceobject

any object with a/b/pa/err_a/err_b properties

AegeanTools.source_finder.get_aux_files(basename)[source]

Look for and return all the aux files that are associated with this filename. Will look for: - background (_bkg.fits) - rms (_rms.fits) - mask (.mim) - catalogue (_comp.fits) - psf map (_psf.fits)

will return filenames if they exist, or None where they do not.

Parameters
basenamestr

The name/path of the input image.

Returns
auxdict

Dict of filenames or None with keys (bkg, rms, mask, cat, psf)

AegeanTools.source_finder.pa_limit(pa)[source]

Position angle is periodic with period 180 deg Constrain pa such that -90<pa<=90

Parameters
pafloat

Initial position angle.

Returns
pafloat

Rotate position angle.

AegeanTools.source_finder.priorized_islands_parinfo(sources, im, wcshelper, stage=3)[source]

Turn a list of sources into a set of islands and parameter estimates which can then be characterised.

Parameters
sources[AegeanTools.models.SimpleSource, … ]

A list of sources in the catalogue.

imnp.ndarray

The image map

wcshelperAegeanTools.wcs_helpers.WCSHelper

A wcs object valid for the image map

stageint

The priorized fitting stage which determines which parameters are fit/fixed. One of: 1 - Fit for flux only. All other params are fixed. 2 - Fit for flux and position. Shape parameters are fixed. 3 - Fit for flux, position, and shape.

Returns
islands[AegeanTools.models.ComponentSource, …]

a list of components

AegeanTools.source_finder.save_catalogue(sources, output, format=None)[source]

Write a catalogue of sources

Parameters
sources[AegeanTools.models.SimpleSource, … ]

A list of characterised sources of type SimpleSource, ComponentSource, or IslandSource.

outputstr

Output filename

formatstr

A descriptor of the output format. Options are: #TODO add a bunch of options ‘auto’ or None - infer from filename extension

Returns
None
AegeanTools.source_finder.theta_limit(theta)[source]

Angle theta is periodic with period pi. Constrain theta such that -pi/2<theta<=pi/2.

Parameters
thetafloat

Input angle.

Returns
thetafloat

Rotate angle.

wcs_helpers

This module contains two classes that provide WCS functions that are not part of the WCS toolkit, as well as some wrappers around the provided tools to make them a lot easier to use.

class AegeanTools.wcs_helpers.Beam(a, b, pa)[source]

Small class to hold the properties of the beam. Properties are a,b,pa. No assumptions are made as to the units, but both a and b have to be >0.

class AegeanTools.wcs_helpers.WCSHelper(wcs, beam, pixscale, refpix, psf_file=None)[source]

A wrapper around astropy.wcs that provides extra functionality, and hides the c/fortran indexing troubles.

Additionally allow psf information to be described in a map instead of the fits header of the image.

Useful functions not provided by astropy.wcs

  • sky2pix/pix2sky functions for vectors and ellipses.

  • functions for calculating the beam in sky/pixel coords

  • the ability to change the beam according to dec-lat

This class tracks both the synthesized beam of the image (beam) and the point spread function (psf). You may think that these things are the same and interchangeable but they are not always. The beam is defined in the wcs of the image header, while the psf can be defined by providing a new image file with 3 dimensions (ra, dec, psf) where the psf is (a, b, pa). # TODO: Check that the above is consistent with the code, and adjust until they are.

Attributes
wcsastropy.wcs.WCS

WCS object

beamAegeanTools.wcs_helpers.Beam

The synthesized beam as defined by the fits header (at the reference location).

pixscale(float, float)

The pixel scale at the reference location (degrees)

refpix(float, float)

The reference location in pixel coordinates

psf_filestr

Filename for a psf map

psf_mapnp.ndarray

A map of the psf as a function of sky position.

psf_wcsnp.ndarray

The WCS object for the psf map

classmethod from_file(filename, beam=None, psf_file=None)[source]

Create a new WCSHelper class from a given fits file.

Parameters
filenamestring

The file to be read

beamAegeanTools.wcs_helpers.Beam or None

The synthesized beam. If the supplied beam is None then one is constructed form the header.

psf_filestr

Filename for a psf map

Returns
objAegeanTools.wcs_helpers.WCSHelper

A helper object

classmethod from_header(header, beam=None, psf_file=None)[source]

Create a new WCSHelper class from the given header.

Parameters
headerastropy.fits.HDUHeader or string

The header to be used to create the WCS helper

beamAegeanTools.wcs_helpers.Beam or None

The synthesized beam. If the supplied beam is None then one is constructed form the header.

psf_filestr

Filename for a psf map

Returns
objAegeanTools.wcs_helpers.WCSHelper

A helper object.

get_beamarea_deg2(ra, dec)[source]

Calculate the area of the synthesized beam in square degrees.

Parameters
ra, decfloat

The sky coordinates at which the calculation is made.

Returns
areafloat

The beam area in square degrees.

get_beamarea_pix(ra, dec)[source]

Calculate the beam area in square pixels.

Parameters
ra, decfloat

The sky coordinates at which the calculation is made

dec
Returns
areafloat

The beam area in square pixels.

get_psf_pix2pix(x, y)[source]

Determine the beam in pixels at the given location in pixel coordinates. The psf is in pixel coordinates.

Parameters
x , yfloat

The pixel coordinates at which the beam is determined.

Returns
a, b, theta(float, float, float)

The psf semi-major axis (pixels), semi-minor axis (pixels), and rotation angle (degrees). If a psf is defined then it is the psf that is returned, otherwise the image restoring beam is returned.

get_psf_sky2pix(ra, dec)[source]

Determine the psf (a,b,pa) at a given sky location. The psf is in pixel coordinates.

Parameters
ra, decfloat

The sky position (degrees).

Returns
a, b, pa(float, float, float)

The psf semi-major axis (pixels), semi-minor axis (pixels), and rotation angle (degrees). If a psf is defined then it is the psf that is returned, otherwise the image restoring beam is returned.

get_psf_sky2sky(ra, dec)[source]

Determine the point spread function in sky coordinates at a given sky location. The psf is returned in degrees.

Parameters
ra, decfloat

The sky position (degrees).

Returns
a, b, pa(float, float, float)

The psf semi-major axis, semi-minor axis, and position angle in (degrees). If a psf is defined then it is the psf that is returned, otherwise the image restoring beam is returned.

get_skybeam(ra, dec)[source]

Determine the beam at the given sky location.

Parameters
ra, decfloat

The sky coordinates at which the beam is determined.

Returns
beamAegeanTools.wcs_helpers.Beam

A beam object, with a/b/pa in sky coordinates

pix2sky(pixel)[source]

Convert pixel coordinates into sky coordinates. Computed on the image wcs.

Parameters
pixel(float, float)

The (x,y) pixel coordinates

Returns
sky(float, float)

The (ra,dec) sky coordinates in degrees

pix2sky_ellipse(pixel, sx, sy, theta)[source]

Convert an ellipse from pixel to sky coordinates.

Parameters
pixel(float, float)

The (x, y) coordinates of the center of the ellipse.

sx, syfloat

The major and minor axes (FHWM) of the ellipse, in pixels.

thetafloat

The rotation angle of the ellipse (degrees). theta = 0 corresponds to the ellipse being aligned with the x-axis.

Returns
ra, decfloat

The (ra, dec) coordinates of the center of the ellipse (degrees).

a, bfloat

The semi-major and semi-minor axis of the ellipse (degrees).

pafloat

The position angle of the ellipse (degrees).

pix2sky_vec(pixel, r, theta)[source]

Given and input position and vector in pixel coordinates, calculate the equivalent position and vector in sky coordinates.

Parameters
pixel(float, float)

origin of vector in pixel coordinates

rfloat

magnitude of vector in pixels

thetafloat

angle of vector in degrees

Returns
ra, decfloat

The (ra, dec) of the origin point (degrees).

r, pafloat

The magnitude and position angle of the vector (degrees).

psf_sky2pix(pos)[source]

Convert sky coordinates into pixel coordinates. Computed on the psf wcs.

Parameters
pos(float, float)

The (ra, dec) sky coordinates (degrees)

Returns
pixel(float, float)

The (x,y) pixel coordinates

sky2pix(pos)[source]

Convert sky coordinates into pixel coordinates. Computed on the image wcs.

Parameters
pos(float, float)

The (ra, dec) sky coordinates (degrees)

Returns
pixel(float, float)

The (x,y) pixel coordinates

sky2pix_ellipse(pos, a, b, pa)[source]

Convert an ellipse from sky to pixel coordinates.

Parameters
pos(float, float)

The (ra, dec) of the ellipse center (degrees).

a, b, pa: float

The semi-major axis, semi-minor axis and position angle of the ellipse (degrees).

Returns
x, yfloat

The (x, y) pixel coordinates of the ellipse center.

sx, syfloat

The major and minor axes (FWHM) in pixels.

thetafloat

The rotation angle of the ellipse (degrees). theta = 0 corresponds to the ellipse being aligned with the x-axis.

sky2pix_vec(pos, r, pa)[source]

Convert a vector from sky to pixel coords. The vector has a magnitude, angle, and an origin on the sky.

Parameters
pos(float, float)

The (ra, dec) of the origin of the vector (degrees).

rfloat

The magnitude or length of the vector (degrees).

pafloat

The position angle of the vector (degrees).

Returns
x, yfloat

The pixel coordinates of the origin.

r, thetafloat

The magnitude (pixels) and angle (degrees) of the vector.

sky_sep(pix1, pix2)[source]

calculate the GCD sky separation (degrees) between two pixels.

Parameters
pix1, pix2(float, float)

The (x,y) pixel coordinates for the two positions.

Returns
distfloat

The distance between the two points (degrees).

AegeanTools.wcs_helpers.fix_aips_header(header)[source]

Search through an image header. If the keywords BMAJ/BMIN/BPA are not set, but there are AIPS history cards, then we can populate the BMAJ/BMIN/BPA. Fix the header if possible, otherwise don’t. Either way, don’t complain.

Parameters
headerastropy.io.fits.HDUHeader

Fits header which may or may not have AIPS history cards.

Returns
headerastropy.io.fits.HDUHeader

A header which has BMAJ, BMIN, and BPA keys, as well as a new HISTORY card.

AegeanTools.wcs_helpers.get_beam(header)[source]

Create a AegeanTools.wcs_helpers.Beam object from a fits header.

BPA may be missing but will be assumed to be zero.

if BMAJ or BMIN are missing then return None instead of a beam object.

Parameters
headerastropy.io.fits.HDUHeader

The fits header.

Returns
beamAegeanTools.wcs_helpers.Beam

Beam object, with a, b, and pa in degrees.

AegeanTools.wcs_helpers.get_pixinfo(header)[source]

Return some pixel information based on the given hdu header pixarea - the area of a single pixel in deg2 pixscale - the side lengths of a pixel (assuming they are square)

Parameters
headerastropy.io.fits.HDUHeader

FITS header information

Returns
pixareafloat

The are of a single pixel at the reference location, in square degrees.

pixscale(float, float)

The pixel scale in degrees, at the reference location.

Notes

The reference location is not always at the image center, and the pixel scale/area may change over the image, depending on the projection.