FilterResponse

class speclite.filters.FilterResponse(wavelength, response, meta, band_shift=None)[source] [edit on github]

Bases: object

A filter response curve tabulated in wavelength.

Some standard filters are included in this package and can be initialized using load_filter(). For example:

>>> rband = load_filter('sdss2010-r')

Objects behave like functions that evaluate their response at aribtrary wavelengths. Wavelength units can be specified, or else default to default_wavelength_unit:

>>> round(rband(6000 * astropy.units.Angstrom), 4)
0.4692
>>> round(rband(6000), 4)
0.4692
>>> round(rband(0.6 * astropy.units.micron), 4)
0.4692

Filters can be also evaluated for an arbitrary array of wavelengths, returning a numpy array of response values:

>>> resp = rband([5980, 6000, 6020])
>>> np.round(resp[1], 4)
0.4692

The effective wavelength of a filter is defined as the photon-weighted mean wavelength:

\[\lambda_{eff} \equiv F[R, \lambda] / F[R, 1]\]

where \(F\) is our convolution operator defined above. Use the effective_wavelength attribute to access this value:

>>> print(np.round(rband.effective_wavelength, 1))
6205.8 Angstrom

The examples below show three different ways to calculate the AB magnitude in the sdss2010-r filter for a source with a constant spectral flux density per unit wavelength of \(10^{-17} \text{erg}/(\text{cm}^2\, \text{s} \,\AA)\). First, we specify the spectrum with a function object:

>>> flux = lambda wlen: 1e-17 * default_flux_unit
>>> print(rband.get_ab_magnitude(flux).round(3))
21.138

Next, we tabulate a constant flux using only two wavelength points that span the filter response:

>>> wlen = [5300, 7200] * default_wavelength_unit
>>> flux = [1e-17, 1e-17] * default_flux_unit
>>> print(rband.get_ab_magnitude(flux, wlen).round(3))
21.138

Since this spectrum undersamples the filter response, it is automatically interpolated. Finally, we tabulate a constant flux using a dense wavelength grid that oversamples the filter response and does not require any interpolation:

>>> wlen = np.linspace(5300, 7200, 200) * default_wavelength_unit
>>> flux = np.ones_like(wlen.value) * 1e-17 * default_flux_unit
>>> print(rband.get_ab_magnitude(flux, wlen).round(3))
21.138

Filters can have an optional wavelength shift applied. See http://arxiv.org/abs/astro-ph/0205243 for applications. Shifted filters can be created from non-shifted filters using create_shifted().

Parameters:
wavelengtharray

A valid array of wavelengths.

responsearray

A dimensionless 1D array of filter response values corresponding to each wavelength. Response values must be non-negative and cannot all be zero. The bounding response values must be zero, and the response is assumed to be zero outside of the specified wavelength range. If this parameter has units, they will be silently ignored.

metadict

A dictionary of metadata which must include values for the keys group_name and band_name, which must be valid python identifiers. However, you are encouraged to provide the full set of keys listed here, and additional keys are also permitted.

band_shiftfloat or None

A shift to apply to filter wavelengths (but not response values). A shifted filter cannot be saved or further shifted, but can otherwise be used like any other filter.

Raises:
ValueError

Invalid wavelength or response input arrays, or missing required keys in the input metadata.

Attributes:
namestr

Canonical name of this filter response in the format “<group_name>-<band_name>”.

effective_wavelengthastropy.units.Quantity

Mean photon-weighted wavelength of this response function, as defined above.

ab_zeropointastropy.units.Quantity

Zeropoint for this filter response in the AB system, as defined above, and including units.

metadict

Dictionary of metadata associated with this filter.

interpolatorscipy.interpolate.interp1d

Linear interpolator of our response function that returns zero for all values outside our wavelength range. Should normally be evaluated through our __call__() convenience method.

Attributes Summary

response

Return the response curve array of the filter

wavelength

Return the wavelength array of the filter

Methods Summary

__call__(wavelength)

Evaluate the filter response at arbitrary wavelengths.

convolve_with_array(wavelength, values[, ...])

Convolve this response with a tabulated function of wavelength.

convolve_with_function(function[, ...])

Convolve this response with a function of wavelength.

create_shifted(band_shift)

Create a copy of this filter response with shifted wavelengths.

get_ab_maggies(spectrum[, wavelength, axis])

Calculate a spectrum's relative AB flux convolution.

get_ab_magnitude(spectrum[, wavelength, axis])

Calculate a spectrum's AB magnitude.

pad_spectrum(spectrum, wavelength[, axis, ...])

Pad a tabulated spectrum to cover this filter's response.

save([directory_name, overwrite])

Save this filter response to file.

Attributes Documentation

response

Return the response curve array of the filter

Returns:
np.ndarray
An array of response values after trimming any extra leading

or trailing zero response values.

wavelength

Return the wavelength array of the filter

Returns:
np.ndarray

An array of wavelengths with their associated units after trimming any extra leading or trailing zero response values.

Methods Documentation

__call__(wavelength)[source] [edit on github]

Evaluate the filter response at arbitrary wavelengths.

Parameters:
wavelengtharray or float

A single wavelength value or an array of wavelengths. If units are included, they will be correctly interpreted. Otherwise default_wavelength_unit is assumed.

Returns:
numpy.ndarray

Numpy array of response values corresponding to each input wavelength.

Raises:
astropy.units.UnitConversionError

Input wavelength(s) have unit that is not convertible to default_wavelength_unit.

convolve_with_array(wavelength, values, photon_weighted=True, interpolate=False, axis=-1, units=None, method='trapz')[source] [edit on github]

Convolve this response with a tabulated function of wavelength.

This is a convenience method that creates a temporary FilterConvolution object to perform the convolution. See that class’ documentation for details on this method’s parameters and usage. See also the notes above about how the convolution integrand is sampled.

Parameters:
wavelengtharray

A valid array of wavelengths that must cover the full range of this filter’s response.

valuesarray or astropy.units.Quantity

Array of function values to use. Values are assumed to be tabulated on our wavelength grid. Values can be multidimensional, in which case an array of convolution results is returned. If values have units, then these are propagated to the result.

photon_weightedbool

Use weights appropriate for a photon-counting detector such as a CCD when this parameter is True. Otherwise, use unit weights.

interpolatebool

Allow interpolation of the tabulated function if necessary. See FilterConvolution for details.

axisint

In case of multidimensional function values, this specifies the index of the axis corresponding to the wavelength dimension.

unitsastropy.units.Quantity or None

When this parameter is not None, then any explicit values units must be convertible to these units, and these units will be applied to the values if they do not already have units.

methodstr

Specifies the numerical integration scheme to use. See FilterConvolution for details.

Returns:
float or array or astropy.units.Quantity

Convolution integral result(s). If the input values have units then these are propagated to the result(s). If the input is multidimensional, then so is the result but with the specified axis integrated out.

convolve_with_function(function, photon_weighted=True, units=None, method='trapz')[source] [edit on github]

Convolve this response with a function of wavelength.

Returns a numerical estimate of the convolution integral \(F[R,f]\) defined above for an arbitrary function of wavelength \(f(\lambda)\). For example, to calculate a filter’s effective wavelength:

>>> rband = load_filter('sdss2010-r')
>>> one = astropy.units.Quantity(1.)
>>> numer = rband.convolve_with_function(lambda wlen: wlen)
>>> denom = rband.convolve_with_function(lambda wlen: one)
>>> print(np.round(numer / denom, 1))
6205.8 Angstrom

Similarly, a filter’s zeropoint can be calculated using:

>>> zpt = rband.convolve_with_function(ab_reference_flux)
>>> print(zpt.round(1))
493486.7 1 / (cm2 s)

Note that both of these values are pre-calculated in the constructor and are available from the effective_wavelength and ab_zeropoint attributes.

Parameters:
functioncallable

Any function that expects a wavelength or array of wavelengths and returns its value. Functions will be called first with wavelength units included and then without units included, in which case they should treat all wavelengths as having default_wavelength_unit. If a function returns a value with units, this will be correctly propagated to the convolution result.

photon_weightedbool

Use weights appropriate for a photon-counting detector such as a CCD when this parameter is True. Otherwise, use unit weights.

unitsastropy.units.Quantity or None

When this parameter is not None, then any explicit units returned by the function must be convertible to these units, and these units will be applied if the function values do not already have units.

methodstr

Specifies the numerical integration scheme to use and must be either ‘trapz’ or ‘simps’, to select the corresponding scipy.integration function. The ‘simps’ method may be more accurate than the default ‘trapz’ method, but should be used with care since it is also less robust and more sensitive to the wavelength grid.

Returns:
float or astropy.units.Quantity

Result of the convolution integral. Units will be included if the function returns a value with units. Otherwise, the implicit units of the result are equal to the implicit function value units multiplied by default_wavelength_unit.

Raises:
ValueError

Function does not behave as expected or invalid method.

RuntimeError

Function returns inconsistent units.

create_shifted(band_shift)[source] [edit on github]

Create a copy of this filter response with shifted wavelengths.

A filter response \(R(\lambda)\) is transformed to blue-shifted response \(R((1+z)\lambda)\) by shifting the wavelengths where its response values are tabulated from \(\lambda_i\) to \(\lambda_i/(1+z)\). Note that only the tabulated wavelengths are transformed and not the response values. See http://arxiv.org/abs/astro-ph/0205243 for applications. A shifted filter cannot be saved or further shifted, but can otherwise be used like any other filter.

Parameters:
band_shiftfloat or None

Shift to apply to filter wavelengths (but not response values).

Returns:
FilterResponse

A new filter object that is a copy of this filter but with the specified wavelength shift applied.

Raises:
RuntimeError

Cannot apply a second wavelength shift to a filter response.

get_ab_maggies(spectrum, wavelength=None, axis=-1)[source] [edit on github]

Calculate a spectrum’s relative AB flux convolution.

The result is the dimensionless ratio \(F[R,f_\lambda] / F[R,f_{\lambda,0}]\) defined above, and provides a linear measure of a source’s flux through this filter relative to an AB standard flux.

Use get_ab_magnitude() for the corresponding AB magnitude.

Parameters:
spectrumcallable or array or astropy.units.Quantity

The spectrum whose flux should be compared with the AB standard flux in this filter band. Can either be a callable object (see convolve_with_function() for details) or else an array (See convolve_with_array() for details). A multidimensional array can be used to calculate maggies for many spectra at once. The spectrum fluxes must either have explicit units that are convertible to default_flux_unit, or else they will be implicitly interpreted as having these default units.

wavelengtharray or astropy.units.Quantity or None

When this parameter is None, the spectrum must be a callable object. Otherwise, the spectrum must be a valid array of wavelengths that must cover the full range of this filter’s response. A spectrum that does not cover this filter can be padded first using pad_spectrum().

axisint

The axis along which wavelength increases in a spectrum array. Ignored when the wavelength parameter is None.

Returns:
float or array
get_ab_magnitude(spectrum, wavelength=None, axis=-1)[source] [edit on github]

Calculate a spectrum’s AB magnitude.

Use get_ab_maggies() for the corresponding dimensionless ratio \(F[R,f_\lambda] / F[R,f_{\lambda,0}]\). The magnitude is calculated as:

Parameters:
spectrumcallable or array or astropy.units.Quantity

See get_ab_maggies() for details.

wavelengtharray or astropy.units.Quantity or None

See get_ab_maggies() for details.

axisint

See get_ab_maggies() for details.

Returns:
float or array
pad_spectrum(spectrum, wavelength, axis=-1, method='median')[source] [edit on github]

Pad a tabulated spectrum to cover this filter’s response.

This is a convenience method that pads an input spectrum (or spectra) so that it covers this filter’s response and can then be used for convolutions and magnitude calculations:

>>> rband = load_filter('decam2014-r')
>>> wlen = np.arange(4000, 10000)
>>> flux = np.ones((100, len(wlen)))
>>> mag = rband.get_ab_magnitude(flux, wlen)
Traceback (most recent call last):
...
ValueError: Wavelengths do not cover decam2014-r
response 3330.0-10990.0 Angstrom.
>>> flux, wlen = rband.pad_spectrum(flux, wlen)
>>> mag = rband.get_ab_magnitude(flux, wlen)

This method does not attempt to perform a physically motivated extrapolation, so should only be used under special circumstances. An appropriate use would be to extend a spectrum to cover wavelengths where the filter response is almost zero or the spectrum is known to be essentially flat.

The three padding methods implemented here are borrowed from numpy.pad():

  • median: Pad with the median flux value.

  • zero: Pad with a zero flux value.

  • edge: Pad with the edge flux values.

If the results from these methods are not in good agreement, you should probably not be padding your spectrum.

Parameters:
spectrumarray or astropy.units.Quantity

See get_ab_maggies() for details. This method does not check for valid flux density units, if any are specified.

wavelengtharray or astropy.units.Quantity or None

See get_ab_maggies() for details.

axisint

See get_ab_maggies() for details.

methodstr

Must be one of ‘median’, ‘zero’, or ‘edge’.

Returns:
tuple

A tuple (padded_spectrum, padded_wavelength) that replaces the inputs with padded equivalents.

save(directory_name='.', overwrite=True)[source] [edit on github]

Save this filter response to file.

The response is saved in the ECSV format and can be read back by passing the returned path to load_filter():

file_name = response.save()
response2 = load_filter(file_name)

The file name in the specified directory will be “<group_name>-<band_name>.ecsv”. Any existing file with the same name will be silently overwritten.

Saved filter responses cannot have any wavelength shift applied since their file naming convention does not distinguish between different shifts. Therefore, attempting to save a filter with a wavelength shift applied will raise a RuntimeError.

Parameters:
directory_namestr

An existing directory where the response file should be written.

overwritebool

Silently overwrite this file when True.

Returns:
str

Absolute path of the created file, including the .ecsv extension.

Raises:
ValueError

Directory name does not exist or refers to a file.

RuntimeError

Filter response has a wavelength shift applied.