Synthetic Photometry (synphot)

Introduction

synphot simulates photometric data and spectra, observed or otherwise. You can incorporate your own filters, spectra, and data. You can also use a pre-defined standard star (Vega), bandpass, or extinction law. Furthermore, it allows you to:

  • Construct complicated composite spectra using different models.
  • Simulate observations.
  • Compute photometric properties such as count rate, effective wavelength, and effective stimulus.
  • Manipulate a spectrum; e.g., applying redshift or normalize it to a given flux value in a given bandpass.
  • Sample a spectrum at given wavelengths.
  • Plot a quick-view of a spectrum.
  • Perform repetitive operations such as simulating the observations of multiple type of sources through multiple bandpasses.

This package covers the general functionalities not related to any particular observatory. If you use HST (and soon JWST), you might also be interested in stsynphot (https://github.com/spacetelescope/stsynphot_refactor), which covers synthetic photometry for the telescope(s).

If you use synphot in your work, please see CITATION for details on how to cite it in your publications.

If you have questions or concerns regarding the software, please open an issue at https://github.com/spacetelescope/synphot_refactor/issues (if not already reported) or contact STScI Help Desk.

Installation and Setup

synphot works under both Python 2 and 3. It requires the following packages:

  • numpy (>= 1.9)
  • astropy (>= 1.3)
  • scipy (>= 0.14)
  • matplotlib (optional for plotting)

You can install synphot using one of the following ways:

  • From AstroConda:

    conda install synphot -c http://ssb.stsci.edu/astroconda
    
  • From standalone release:

    pip install git+https://github.com/spacetelescope/synphot_refactor.git@0.1
    
  • From source (example shown is for the dev version):

    git clone https://github.com/spacetelescope/synphot_refactor.git
    cd synphot_refactor
    python setup.py install
    

To use the pre-defined standard star, extinction laws, and bandpasses, it is recommended that you copy synphot.cfg to your $HOME/.astropy/config/ directory. In doing so, you can avoid connecting directly to STScI HTTP service, which is slower and might not be available all the time.

If necessary (i.e., if you are not an internal STScI user), replace /grp/hst/cdbs with your local directory name (e.g., /my/local/dir/cdbs), where you plan to store the data files. Then, follow the instructions below to “install” the data files.

Note

synphot data files are a minimal subset of those required by stsynphot. If you plan to use the latter anyway, please also read the instructions in its documentation.

To install local data files via HTTP into, say, /my/local/dir/cdbs:

>>> from synphot.util import download_data
>>> file_list = download_data('/my/local/dir/cdbs')

If you have your own version of the data files that you wish to use, you can modify synphot.cfg to point to your own copies without having to download using the function above. However, before you do so, make sure that your own file(s) can be read in successfully with read_fits_spec(). For example, if you want to use your own Johnson V throughput file, you can modify this line in your $HOME/.astropy/config/synphot.cfg file:

johnson_v_file = /my/other/dir/my_johnson_v.fits

Alternately, you can also take advantage of Configuration system (astropy.config) to manage synphot data files. For example, you can also temporarily use a different Johnson V throughput file:

>>> from synphot.config import conf
>>> print(conf.johnson_v_file)
/my/local/dir/cdbs/comp/nonhst/johnson_v_004_syn.fits
>>> with conf.set_temp('johnson_v_file', '/my/other/dir/my_johnson_v.fits'):
...     print(conf.johnson_v_file)
/my/other/dir/my_johnson_v.fits
>>> print(conf.johnson_v_file)
/my/local/dir/cdbs/comp/nonhst/johnson_v_004_syn.fits

Getting Started

This section only contains minimal examples showing how to use this package. For detailed documentation, see Using synphot.

In the examples below, you will notice that most models are from synphot.models, not astropy.modeling.models, because the models in synphot have extra things like sampleset that are not (yet) available in Astropy. Despite this, some models like Const1D does not need the extra things to work, so they can be used directly. When in doubt, see if a model is in synphot.models first before using Astropy’s.

from astropy import units as u
from synphot import units, SourceSpectrum
from synphot.models import (BlackBodyNorm1D, GaussianAbsorption1D,
                            GaussianFlux1D)
from synphot.spectrum import BaseUnitlessSpectrum
# Create a Gaussian absorption profile that absorps 80% of the flux at
# 4000 Angstrom with a sigma of 20 Angstrom
g_abs = BaseUnitlessSpectrum(GaussianAbsorption1D, amplitude=0.8,
                             mean=4000, stddev=20)
# Create a Gaussian emission line with the given total flux
# centered at 3000 Angstrom with FWHM of 100 Angstrom
g_em = SourceSpectrum(GaussianFlux1D,
                      total_flux=3.5e-13*u.erg/(u.cm**2 * u.s),
                      mean=3000, fwhm=100)
# Create a blackbody source spectrum with a temperature of 6000 K
bb = SourceSpectrum(BlackBodyNorm1D, temperature=6000)
# Combine the above components to create a source spectrum that is twice
# the original blackbody flux with the Gaussian emission and absorption
# lines
sp = g_abs * (g_em + 2 * bb)
# Plot the spectrum, zooming in on the line features
sp.plot(left=1, right=7000)

(Source code)

_images/index-1.png

Sample the spectrum at 0.3 micron:

>>> sp(0.3 * u.micron)
<Quantity 0.0012685 PHOTLAM>

Or sample the same thing but in a different flux unit:

>>> sp(0.3 * u.micron, flux_unit=units.FLAM)
<Quantity 8.57716634e-15 FLAM>

Sample the spectrum at its native wavelength set:

>>> sp(sp.waveset)
<Quantity [0.00000000e+00, 0.00000000e+00, 0.00000000e+00, ...,
           4.47483287e-05, 4.34264666e-05, 4.21423394e-05] PHOTLAM>

Models that built the spectrum:

>>> print(sp)
SourceSpectrum at z=0.0
Model: CompoundModel2
Inputs: ('x',)
Outputs: ('y',)
Model set size: 1
Expression: ([0] + ([1] | [2])) * [3]
Components:
    [0]: <GaussianFlux1D(amplitude=...>

    [1]: <BlackBodyNorm1D(temperature=6000.0)>

    [2]: <Scale(factor=2.0)>

    [3]: <GaussianAbsorption1D(amplitude=0.8, mean=4000.0, stddev=20.0)>
Parameters:
    ...

Redshift the source spectrum by \(z = 0.2\):

>>> sp.z = 0.2

Create a box-shaped bandpass centered at 4000 Angstrom with a width of 2000 Angstrom:

>>> from synphot import SpectralElement
>>> from synphot.models import Box1D
>>> bp = SpectralElement(Box1D, amplitude=1, x_0=4000, width=2000)

Normalize the source spectrum to 1 Jy in a given box bandpass and integrate it:

>>> sp_rn = sp.normalize(1 * u.Jy, band=bp)
>>> sp_rn.integrate()
<Quantity 12901.13466265 ph / (cm2 s)>

Create an observation by passing the redshifted and normalized source spectrum through the box bandpass:

>>> from synphot import Observation
>>> obs = Observation(sp_rn, bp)

Calculate the count rate of the observation above for an 2-meter telescope:

>>> import numpy as np
>>> area = np.pi * (1 * u.m) ** 2
>>> area
<Quantity 3.14159265 m2>
>>> obs.countrate(area=area)
<Quantity 24219651.1854991 ct / s>

A Brief History

A brief history: First, there was STSDAS SYNPHOT. Then, there was ASTROLIB PYSYNPHOT (Lim et al. 2015), which aimed at replacing STSDAS SYNPHOT using Python. In order to take advantage of Models and Fitting (astropy.modeling) and Units and Quantities (astropy.units) and to repurpose the functionality for a wider audience (other than HST users), it was refactored again and separated into synphot and stsynphot (see Synthetic Photometry (synphot)).

API

synphot.binning Module

Utilities related to wavelength bin calculations.

Functions

calculate_bin_edges(centers) Calculate the edges of wavelength bins given the centers.
calculate_bin_widths(edges) Calculate the widths of wavelengths bins given their edges.
calculate_bin_centers(edges) Calculate the centers of wavelengths bins given their edges.
wave_range(bins, cenwave, npix[, mode]) Calculate the wavelength range covered by the given number of pixels centered on the given central wavelength of the given bins.
pixel_range(bins, waverange[, mode]) Calculate the number of pixels within the given wavelength range and the given bins.

Also imports this C-extension to local namespace:

synphot.config Module

Synphot configurable items.

The default configuration heavily depends on STScI CDBS structure but it can be easily re-configured as the user wishes via astropy.config.

synphot.exceptions Module

Exceptions specific to synthetic photometry.

Classes

SynphotError Base class for synphot exceptions.
TableFormatError(msg[, rows]) Exceptions to do with table access.
DuplicateWavelength(msg[, rows]) Duplicate wavelengths are not allowed in table.
ZeroWavelength(msg[, rows]) Zero wavelengths are not allowed in table.
UnsortedWavelength(msg[, rows]) Unsorted wavelengths are not allowed in table.
OverlapError Exceptions to do with overlap checking.
PartialOverlap Partial overlap is not allowed.
DisjointError Disjoint data is not allowed.
UndefinedBinset Exceptions for undefined bin set.
IncompatibleSources Two sources in composite spectrum are not compatible.
InterpolationNotAllowed Exceptions for interpolation.
ExtrapolationNotAllowed Exceptions for extrapolation.

Class Inheritance Diagram

Inheritance diagram of synphot.exceptions.SynphotError, synphot.exceptions.TableFormatError, synphot.exceptions.DuplicateWavelength, synphot.exceptions.ZeroWavelength, synphot.exceptions.UnsortedWavelength, synphot.exceptions.OverlapError, synphot.exceptions.PartialOverlap, synphot.exceptions.DisjointError, synphot.exceptions.UndefinedBinset, synphot.exceptions.IncompatibleSources, synphot.exceptions.InterpolationNotAllowed, synphot.exceptions.ExtrapolationNotAllowed

synphot.models Module

Spectrum models not in astropy.modeling.

Functions

get_waveset(model) Get optimal wavelengths for sampling a given model.
get_metadata(model) Get metadata for a given model.

Classes

BlackBody1D(*args, **kwargs) Create a blackbody spectrum model with given temperature.
BlackBodyNorm1D(*args, **kwargs) Create a normalized blackbody spectrum with given temperature.
Box1D([amplitude, x_0, width]) Same as astropy.modeling.functional_models.Box1D, except with sampleset defined.
ConstFlux1D(amplitude, **kwargs) One dimensional constant flux model.
Empirical1D(**kwargs) Empirical (sampled) spectrum or bandpass model.
Gaussian1D([amplitude, mean, stddev]) Same as astropy.modeling.functional_models.Gaussian1D, except with sampleset defined.
GaussianAbsorption1D([amplitude, mean, stddev]) Same as astropy.modeling.functional_models.GaussianAbsorption1D, except with sampleset defined.
GaussianFlux1D(*args, **kwargs) Same as Gaussian1D but accepts extra keywords below.
Lorentz1D([amplitude, x_0, fwhm]) Same as astropy.modeling.functional_models.Lorentz1D, except with sampleset defined.
MexicanHat1D([amplitude, x_0, sigma]) Same as astropy.modeling.models.MexicanHat1D, except with sampleset defined.
PowerLawFlux1D(amplitude, x_0, alpha, **kwargs) One dimensional power law model with proper flux handling.
Trapezoid1D([amplitude, x_0, width, slope]) Same as astropy.modeling.functional_models.Trapezoid1D, except with sampleset defined.

Class Inheritance Diagram

Inheritance diagram of synphot.models.BlackBody1D, synphot.models.BlackBodyNorm1D, synphot.models.Box1D, synphot.models.ConstFlux1D, synphot.models.Empirical1D, synphot.models.Gaussian1D, synphot.models.GaussianAbsorption1D, synphot.models.GaussianFlux1D, synphot.models.Lorentz1D, synphot.models.MexicanHat1D, synphot.models.PowerLawFlux1D, synphot.models.Trapezoid1D

synphot.observation Module

This module defines an observed spectrum, i.e., a source spectrum that has gone through a bandpass.

Classes

Observation(spec, band[, binset, force]) This is an observed spectrum, where a source spectrum has gone through a bandpass.

Class Inheritance Diagram

Inheritance diagram of synphot.observation.Observation

synphot.reddening Module

This module defines reddening laws and extinction curves.

Functions

etau_madau(wave, z, **kwargs) Madau 1995 extinction for a galaxy at given redshift.

Classes

ExtinctionModel1D(**kwargs) Model to handle extinction curve.
ReddeningLaw(modelclass[, clean_meta]) Class to handle reddening law.
ExtinctionCurve(modelclass[, clean_meta]) Class to handle extinction curve.

Class Inheritance Diagram

Inheritance diagram of synphot.reddening.ExtinctionModel1D, synphot.reddening.ReddeningLaw, synphot.reddening.ExtinctionCurve

synphot.specio Module

This modules handles synthetic photometry data formats.

Functions

read_remote_spec(filename[, encoding, …]) Read FITS or ASCII spectrum from a remote location.
read_spec(filename[, fname]) Read FITS or ASCII spectrum.
read_ascii_spec(filename[, wave_unit, flux_unit]) Read ASCII spectrum.
read_fits_spec(filename[, ext, wave_col, …]) Read FITS spectrum.
write_fits_spec(filename, wavelengths, fluxes) Write FITS spectrum.

synphot.spectrum Module

This module defines the different types of spectra.

Classes

BaseSpectrum(modelclass[, clean_meta]) Base class to handle spectrum or bandpass.
BaseSourceSpectrum(modelclass[, clean_meta]) Base class to handle spectrum with flux unit like source spectrum and observation.
SourceSpectrum(modelclass[, z, z_type]) Class to handle source spectrum.
BaseUnitlessSpectrum(modelclass[, clean_meta]) Base class to handle unitless spectrum like bandpass, reddening, etc.
SpectralElement(modelclass[, clean_meta]) Class to handle instrument filter bandpass.

synphot.thermal Module

This module defines thermal spectra.

Classes

ThermalSpectralElement(modelclass, temperature) Class to handle spectral element with associated thermal properties.

Class Inheritance Diagram

Inheritance diagram of synphot.thermal.ThermalSpectralElement

synphot.units Module

This module handles photometry units that are not in astropy.units.

Functions

spectral_density_vega(wav, vegaflux) Flux equivalencies between PHOTLAM and VEGAMAG.
spectral_density_count(wav, area) Flux equivalencies between PHOTLAM and count/OBMAG.
convert_flux(wavelengths, fluxes, …) Perform conversion for supported flux units.
validate_unit(input_unit) Validate unit.
validate_wave_unit(wave_unit) Like validate_unit() but specific to wavelength.
validate_quantity(input_value, output_unit) Validate quantity (value and unit).

synphot.utils Module

Synthetic photometry utility functions.

Functions

overlap_status(a, b) Check overlap between two arrays.
validate_totalflux(totalflux) Check integrated flux for invalid values.
validate_wavelengths(wavelengths) Check wavelengths for synphot compatibility.
generate_wavelengths([minwave, maxwave, …]) Generate wavelength array to be used for spectrum sampling.
merge_wavelengths(waveset1, waveset2[, …]) Return the union of the two sets of wavelengths using numpy.union1d().
download_data(cdbs_root[, verbose, dry_run]) Download CDBS data files to given root directory.

References

Calzetti, D., Armus, L., Bohlin, R. C., Kinney, A. L., Koornneef, J., & Storchi-Bergmann, T. 2000, ApJ, 533, 682

Cardelli, J. A., Clayton, G. C., & Mathis, J. S. 1989, ApJ, 345, 245

Gordon, K. D., Clayton, G. C., Misselt, K. A., Landolt, A. U., & Wolff, M. J. 2003, ApJ, 594, 279

Horne, K. 1988, in New Directions in Spectophotometry: A Meeting Held in Las Vegas, NV, March 28-30, Application of Synthetic Photometry Techniques to Space Telescope Calibration, ed. A. G. Davis Philip, D. S. Hayes, & S. J. Adelman (Schenectady, NY: L. Davis Press), 145

Koornneef, J., Bohlin, R., Buser, R., Horne, K., & Turnshek, D. 1986, Highlights Astron., 7, 833

Laidler, V., et al. 2008, Synphot Data User’s Guide, Version 1.2 (Baltimore, MD: STScI)

Lim, P. L., Diaz, R. I., & Laidler, V. 2015, PySynphot User’s Guide (Baltimore, MD: STScI), https://pysynphot.readthedocs.io/en/latest/

Madau, P., et al. 1995, ApJ, 441, 18

Oke, J. B., 1974, ApJS, 27, 21

Rybicki, G. B., & Lightman, A. P. 1979, Radiative Processes in Astrophysics (New York, NY: Wiley)

Schneider, D. P., Gunn, J. E., & Hoessel J. G. 1983, ApJ, 264, 337