Measurements: How you put observations to the Toolbox

To use PDR Toolbox, you need to create Measurements from your observations. A Measurement consists of a value and an error. These can be single-valued or an array of values. In the typical case of an image, the Measurement is a representation of a FITS file with two HDUs, the first HDU is the spatial map of intensity and the 2nd HDU is the spatial map of the errors. It is based on astropy’s CCDData if you are familiar with that. Typical sub-millimeter maps we get from telescopes don’t have the error plane, but PDRT makes it easy for you to create one if you know the magnitude of the error. Typical FITS images will be in intensity units, equivalent to \({\rm erg~s^{-1}~cm^{-2}~sr^{-1}}\), or in \({\rm K~km~s^{-1}}\). For the latter, PDRT will do appropriate conversion as necessary when it uses your images (the original Measurement remains unchanged).

For example how to use Measurements, see the notebook PDRT_Example_Measurements.ipynb.


Manage spectral line or continuum observations

class pdrtpy.measurement.Measurement(*args, **kwargs)[source]

Bases: CCDData

Measurement represents one or more observations of a given spectral line or continuum. It is made up of a value array, an uncertainty array, units, and a string identifier. It is based on astropy.nddata.CCDData. It can represent a single pixel observation or an image. Mathematical operations using Measurements will correctly propagate errors.

Typically, Measurements will be instantiated from a FITS file by using the the read() or make_measurement() methods. For a list of recognized spectral line identifiers, see supported_lines().

Parameters:
Raises:

TypeError – if beam parameters are not Quantities

Measurements can also be instantiated by the read(\*args, \**kwargs), to create an Measurement instance based on a FITS file. This method uses fits_measurement_reader() with the provided parameters. Example usage:

from pdrtpy.measurement import Measurement

my_obs = Measurement.read("file.fits",identifier="CII_158")
my_other_obs = Measurement.read("file2.fits",identifier="CO2_1",
                                 unit="K km/s",
                                 bmaj=9.3*u.arcsec,
                                 bmin=14.1*u.arcsec,
                                 bpa=23.2*u.degrees)

By default image axes with only a single dimension are removed on read. If you do not want this behavior, used read(squeeze=False). See also: astropy.nddata.CCDData.

Attributes:
SN

Return the signal to noise ratio (value/error)

beam

Return the beam parameters as astropy Quantities or None if beam is not set

data

~numpy.ndarray-like : The stored dataset.

dtype

numpy.dtype of this object’s data.

error

Return the underlying error array

filename

The FITS file that created this measurement, or None if it didn’t originate from a file

flags
header
id

Return the string ID of this measurement, e.g., CO_10

levels
mask

any type : Mask for the dataset, if any.

meta

dict-like : Additional meta information about the dataset.

ndim

integer dimensions of this object’s data

psf

Image representation of the PSF for the dataset.

shape

shape tuple of this object’s data.

size

integer size of this object’s data.

title

A formatted title (e.g., LaTeX) that can be in plotting.

uncertainty

any type : Uncertainty in the dataset, if any.

unit

~astropy.units.Unit : Unit for the dataset, if any.

value

Return the underlying data array

wcs

any type : A world coordinate system (WCS) for the dataset, if any.

Methods

add(other)

Add this Measurement to another, propagating errors, units, and updating identifiers.

convert_unit_to(unit[, equivalencies])

Returns a new NDData object whose values have been converted to a new unit.

copy()

Return a copy of the CCDData object.

divide(other)

Divide this Measurement by another, propagating errors, units, and updating identifiers.

from_table(filename[, format, array])

Table file reader for Measurement class.

get(world_x, world_y[, log])

Get the value(s) at the give world coordinates

get_pixel(world_x, world_y)

Return the nearest pixel coordinates to the input world coordinates

identifier(id)

Set the string ID of this measurement, e.g., CO_10

is_ratio()

Indicate if this Measurement is a ratio.

is_single_pixel()

Is this Measurement a single value? :returns: True if a single value (pixel) :rtype: bool

make_measurement(datafile, error, outfile[, ...])

Create a FITS files with 2 HDUS, the first being the datavalue and the 2nd being the data uncertainty.

multiply(other)

Multiply this Measurement by another, propagating errors, units, and updating identifiers.

read

subtract(other)

Subtract another Measurement from this one, propagating errors, units, and updating identifiers.

to_hdu([hdu_mask, hdu_uncertainty, ...])

Creates an HDUList object from a CCDData object.

write(filename, **kwd)

Write this Measurement to a FITS file with value in 1st HDU and error in 2nd HDU.

property SN

Return the signal to noise ratio (value/error)

Return type:

numpy.ndarray

add(other)[source]

Add this Measurement to another, propagating errors, units, and updating identifiers. Masks are logically or’d.

Parameters:

other (Measurement or number) – a Measurement or number to add

property beam

Return the beam parameters as astropy Quantities or None if beam is not set

divide(other)[source]

Divide this Measurement by another, propagating errors, units, and updating identifiers. Masks are logically or’d.

Parameters:

other (Measurement or number) – a Measurement or number to divide by

property error

Return the underlying error array

Return type:

numpy.ndarray

property filename

The FITS file that created this measurement, or None if it didn’t originate from a file

Return type:

str or None

static from_table(filename, format='ipac', array=False)[source]

Table file reader for Measurement class. Create one or more Measurements from a table. The input table header must contain the columns:

data - the data value

uncertainty - the error on the data, can be absolute error or percent. If percent, the header unit row entry for this column must be “%”

identifier - the identifier of this Measurement which should match a model in the ModelSet you are using, e.g., “CII_158” for [C II] 158 $\mu$m

The following columns are optional:

bmaj - beam major axis size

bmin - beam minor axis size

bpa - beam position angle

The table must specify the units of each column, e.g. a unit row in the header for IPAC format. Leave column entry blank if unitless. Units of value and error should be the same or conformable. Units must be transformable to a valid astropy.unit.Unit.

Parameters:
  • filename (str) – Name of table file.

  • formatAstropy Table format format. e.g., ascii, ipac, votable. Default is IPAC format

  • array (bool) – Controls whether a list of Measurements or a single Measurement is returned. If array is True, one Measurement instance will be created for each row in the table and a Python list of Measurements will be returned. If array is False, one Measurement containing all the points in the data member will be returned. If array is False, the identifier and beam parameters of the first row will be used. If feeding the return value to a plot method such as phasespace(), choose array=False. Default:False.

Return type:

Measurement or list of Measurement

get(world_x, world_y, log=False)[source]

Get the value(s) at the give world coordinates

Parameters:
  • world_x (float or array-like) – the x value in world units of naxis1

  • world_y (float or array-lke) – the y value in world units of naxis2

  • log (bool) – True if the input coords are logarithmic Default:False

Returns:

The value(s) of the Measurement at input coordinates

Return type:

float

get_pixel(world_x, world_y)[source]

Return the nearest pixel coordinates to the input world coordinates

Parameters:
  • world_x (float) – The horizontal world coordinate

  • world_y (float) – The vertical world coordinate

property id

Return the string ID of this measurement, e.g., CO_10

Return type:

str

identifier(id)[source]

Set the string ID of this measurement, e.g., CO_10

Parameters:

id (str) – the identifier

is_ratio()[source]

Indicate if this Measurement is a ratio.. This method looks for the ‘/’ past the first character of the` Measurement` identifier, such as “CII_158/CO_32” See also pdrutils.is_ratio(string)

Returns:

True if the Measurement is a ratio, False otherwise

Return type:

bool

is_single_pixel()[source]

Is this Measurement a single value? :returns: True if a single value (pixel) :rtype: bool

property levels
static make_measurement(datafile, error, outfile, rms=None, masknan=True, overwrite=False, unit='adu')[source]

Create a FITS files with 2 HDUS, the first being the datavalue and the 2nd being the data uncertainty. This format makes allows the resulting file to be read into the underlying :class:’~astropy.nddata.CCDData` class.

Parameters:
  • datafile (str) – The FITS file containing the data as a function of spatial coordinates

  • error (str) –

    The errors on the data Possible values for error are:

    • a filename with the same shape as datafile containing the error values per pixel

    • a percentage value ‘XX%’ must have the “%” symbol in it

    • ’rms’ meaning use the rms parameter if given, otherwise look for the RMS keyword in the FITS header of the datafile

  • outfile (str) – The output file to write the result in (FITS format)

  • rms (float or astropy.units.Unit) – If error == ‘rms’, this value may give the rms in same units as data (e.g ‘erg s-1 cm-2 sr-1’).

  • masknan (bool) – Whether to mask any pixel where the data or the error is NaN. Default:true

  • overwrite (bool) – If True, overwrite the output file if it exists. Default: False.

  • unit (astropy.units.Unit or str) – Intensity unit to use for the data, this will override BUNIT in header if present.

Raises:
  • Exception – on various FITS header issues

  • OSError – if overwrite is False and the output file exists.

Example usage:

# example with percentage error
Measurement.make_measurement("my_infile.fits",error='10%',outfile="my_outfile.fits")

# example with measurement in units of K km/s and error
# indicated by RMS keyword in input file.
Measurement.make_measurement("my_infile.fits",error='rms',outfile="my_outfile.fits",unit="K km/s",overwrite=True)
multiply(other)[source]

Multiply this Measurement by another, propagating errors, units, and updating identifiers. Masks are logically or’d.

Parameters:

other (Measurement or number) – a Measurement or number to multiply

subtract(other)[source]

Subtract another Measurement from this one, propagating errors, units, and updating identifiers. Masks are logically or’d.

Parameters:

other (Measurement or number) – a Measurement or number to subtract

property title

A formatted title (e.g., LaTeX) that can be in plotting.

Return type:

str or None

property value

Return the underlying data array

Return type:

numpy.ndarray

write(filename, **kwd)[source]

Write this Measurement to a FITS file with value in 1st HDU and error in 2nd HDU. See astropy.nddata.CCDData.write().

Parameters:
  • filename (str) – Name of file.

  • kwd – All additional keywords are passed to astropy.io.fits

The available built-in formats are:

Format

Read

Write

Auto-identify

pdrtpy.measurement.fits_measurement_reader(filename, hdu=0, unit=None, hdu_mask='MASK', hdu_flags=None, key_uncertainty_type='UTYPE', **kwd)[source]

FITS file reader for Measurement class, which will be called by Measurement.read().

Parameters:
  • filename (str) – Name of FITS file.

  • identifier (str) – string indicating what this is an observation of, e.g., “CO_10” for CO(1-0)

  • squeeze (bool) – If True, remove single dimension axes from the input image. Default: True

  • hdu (int, optional) – FITS extension from which Measurement should be initialized. If zero and and no data in the primary extension, it will search for the first extension with data. The header will be added to the primary header. Default is 0.

  • unit (astropy.units.Unit, optional) – Units of the image data. If this argument is provided and there is a unit for the image in the FITS header (the keyword BUNIT is used as the unit, if present), this argument is used for the unit. Default is None.

  • hdu_uncertainty (str or None, optional) – FITS extension from which the uncertainty should be initialized. If the extension does not exist the uncertainty of the Measurement is None. Default is 'UNCERT'.

  • hdu_mask (str or None, optional) – FITS extension from which the mask should be initialized. If the extension does not exist the mask of the Measurement is None. Default is 'MASK'.

  • hdu_flags (str or None, optional) – Currently not implemented. Default is None.

  • kwd – Any additional keyword parameters are passed through to the FITS reader in astropy.io.fits

Raises:

TypeError – If the conversion from CCDData to Measurement fails