Reference/API
pysynphot.binning
Utilities related to wavelength bin calculations.
- pysynphot.binning.calculate_bin_centers(edges)
Calculate the centers of wavelengths bins given their edges.
- Parameters:
edges (array_like) – Sequence of bin edges. Must be 1D and have at least two values.
- Returns:
centers – Array of bin centers. Will be 1D and have one less value than
edges
.- Return type:
ndarray
- pysynphot.binning.calculate_bin_edges(centers)
Calculate the edges of wavelength bins given the centers.
The algorithm calculates bin edges as the midpoints between bin centers and treats the first and last bins as symmetric about their centers.
- Parameters:
centers (array_like) – Sequence of bin centers. Must be 1D and have at least two values.
- Returns:
edges – Array of bin edges. Will be 1D and have one more value than
centers
.- Return type:
ndarray
- pysynphot.binning.calculate_bin_widths(edges)
Calculate the widths of wavelengths bins given their edges.
- Parameters:
edges (array_like) – Sequence of bin edges. Must be 1D and have at least two values.
- Returns:
widths – Array of bin widths. Will be 1D and have one less value than
edges
.- Return type:
ndarray
pysynphot.Cache
This module is a container for IO-intensive items that should be read in only once, and then re-used from memory.
This includes the reddening laws
(pysynphot.locations.RedLaws
)
and some indices for the catalog
model atlases
(pysynphot.Cache.CATALOG_CACHE
).
- pysynphot.Cache.reset_catalog_cache()
Empty the
CATALOG_CACHE
global variable.
pysynphot.catalog
This module is useful for working with catalog spectra such as Castelli-Kurucz Atlas, Kurucz Atlas, and Phoenix Models.
Spectra are constructed from basis spectra which are indexed for various
combinations of effective temperature (\(T_{\mathrm{eff}}\)),
metallicity ([M/H]
), and log surface gravity (\(\log g\)).
The user may specify any combination of \(T_{\mathrm{eff}}\),
[M/H]
, and \(\log g\) so long as each parameter is within the range
for that parameter defined by the catalog.
For example, the Castelli-Kurucz Atlas catalog contains spectra for effective temperatures between 3500 and 50000 K. In this case, no spectrum can be constructed for \(T_{\mathrm{eff}}\) of 50001 K or 3499 K. The range of parameters available may vary from catalog to catalog.
More information on catalogs can be found in Appendix A: Catalogs and Spectral Atlases.
- class pysynphot.catalog.Icat(catdir, Teff, metallicity, log_g)
Bases:
TabularSourceSpectrum
This class constructs a model from the grid available in catalogs. Specifically, they are Castelli-Kurucz Atlas, Kurucz Atlas, and Phoenix Models.
Each grid contains a master file named “catalog.fits”, as defined by
pysynphot.locations.CAT_TEMPLATE
. The basis spectra are located atpysynphot.locations.KUR_TEMPLATE
. You may inspect the data files in CRDS to see how they are formatted.- Parameters:
- waveunits, fluxunits
Catalog units for wavelength and flux.
- Type:
- wave, flux
Wavelength set and associated flux in catalog units.
- Type:
array_like
- Raises:
pysynphot.exceptions.ParameterOutOfBounds – Given parameter value is out of bounds.
Examples
>>> spec = S.Icat('k93models', 6440, 0, 4.3)
pysynphot.exceptions
Custom exceptions for pysynphot
to raise.
- class pysynphot.exceptions.PysynphotError(msg)
Bases:
Exception
Parent class for
pysynphot
exceptions.- Parameters:
msg (str) – Error message.
- class pysynphot.exceptions.TableFormatError(msg, rows=None)
Bases:
PysynphotError
Exception to do with table access.
- class pysynphot.exceptions.DuplicateWavelength(msg, rows=None)
Bases:
TableFormatError
Exception for duplicate wavelength values in table.
- class pysynphot.exceptions.ZeroWavelength(msg, rows=None)
Bases:
TableFormatError
Exception for wavelength values containing zero.
- class pysynphot.exceptions.UnsortedWavelength(msg, rows=None)
Bases:
TableFormatError
Exception for wavelength values not in ascending or descending order.
- class pysynphot.exceptions.BadRow(msg, rows=None)
Bases:
TableFormatError
Exception for invalid row in table.
- class pysynphot.exceptions.OverlapError(msg)
Bases:
PysynphotError
Exception to do with overlap checking.
- class pysynphot.exceptions.PartialOverlap(msg)
Bases:
OverlapError
Exception for partial overlap between two spectra.
- class pysynphot.exceptions.DisjointError(msg)
Bases:
OverlapError
Exception for no overlap between two spectra.
- class pysynphot.exceptions.GraphtabError(msg)
Bases:
PysynphotError
Exception to do with graph table traversal.
- class pysynphot.exceptions.UnusedKeyword(msg)
Bases:
GraphtabError
Exception for unused keyword in graph table lookup.
- class pysynphot.exceptions.IncompleteObsmode(msg)
Bases:
GraphtabError
Exception for incomplete observation mode in graph table lookup.
- class pysynphot.exceptions.AmbiguousObsmode(msg)
Bases:
GraphtabError
Exception for ambiguous observation mode in graph table lookup.
- class pysynphot.exceptions.UndefinedBinset(msg)
Bases:
PysynphotError
Exception for undefined
binset
in bandpass or observation.
- class pysynphot.exceptions.ExtrapolationNotAllowed(msg)
Bases:
PysynphotError
Exception for invalid extrapolation.
- class pysynphot.exceptions.ParameterOutOfBounds(msg)
Bases:
PysynphotError
Exception for invalid parameter value in a catalog.
- class pysynphot.exceptions.IncompatibleSources(msg)
Bases:
PysynphotError
Exception for operation on two incompatible spectra types.
pysynphot.extinction
This module handles deprecated extinction models for backward compatibility with IRAF STSDAS SYNPHOT.
- class pysynphot.extinction.DeprecatedExtinction(extval, redlaw)
Bases:
SpectralElement
This class handles deprecated extinction models from IRAF STSDAS SYNPHOT like a spectral element.
- Parameters:
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- wave, throughput
Wavelength set in
waveunits
and associated unitless extinction.- Type:
array_like
Examples
>>> extinction = S.Extinction(0.3, 'gal1')
- class pysynphot.extinction.Gal1(extval)
Deprecated Milky Way extinction curve (Seaton 1979).
- Parameters:
extval (float) – Value of \(E(B-V)\) in magnitudes.
- transparencytable
This is the same as \(\mathrm{THRU}\) defined in
reddening()
.- Type:
array_like
- class pysynphot.extinction.Smc(extval)
Deprecated SMC extinction curve (Prevot et al. 1984).
- Parameters:
extval (float) – Value of \(E(B-V)\) in magnitudes.
- transparencytable
This is the same as \(\mathrm{THRU}\) defined in
reddening()
.- Type:
array_like
- class pysynphot.extinction.Lmc(extval)
Deprecated LMC extinction curve (Howarth 1983).
- Parameters:
extval (float) – Value of \(E(B-V)\) in magnitudes.
- transparencytable
This is the same as \(\mathrm{THRU}\) defined in
reddening()
.- Type:
array_like
- class pysynphot.extinction.Xgal(extval)
Deprecated Extra-galactic extinction curve (Calzetti et al. 1994).
- Parameters:
extval (float) – Value of \(E(B-V)\) in magnitudes.
- transparencytable
This is the same as \(\mathrm{THRU}\) defined in
reddening()
.- Type:
array_like
pysynphot.locations
This module handles locations of data files.
Global Variables
pysynphot.locations.rootdir
- Root directory for TRDS/CRDS data files. By default, it is extracted from yourPYSYN_CDBS
environment variable.pysynphot.locations.specdir
- Data directory for data files distributed with this software.pysynphot.locations.CAT_TEMPLATE
andpysynphot.locations.KUR_TEMPLATE
- String templates used forIcat
to select spectra from catalogs.pysynphot.locations.VegaFile
- Vega spectrum to use forvegamag
calculations.pysynphot.locations.EXTDIR
- Directory containing extinction curves.pysynphot.locations.RedLaws
- Dictionary mapping reddening laws to data files or cached instances (seeCache
).pysynphot.locations.wavecat
- Data file forwavetable
.pysynphot.locations.CONVERTDICT
- Dictionary mapping IRAF-style directory shortcuts to actual paths.
- pysynphot.locations.get_data_filename(filename)
Map filename to its actual path.
- pysynphot.locations.get_latest_file(template, raise_error=False, err_msg='')
Find the filename that appears last in sorted order based on given template.
- Parameters:
template (str) – Search template in the form of
path/pattern
where pattern is acceptable byfnmatch
.raise_error (bool, optional) – Raise an error when no files found. Otherwise, will issue warning only.
err_msg (str) – Alternate message for when no files found. If not given, generic message is used.
- Returns:
filename – Latest filename.
- Return type:
- Raises:
IOError – No files found.
- pysynphot.locations.irafconvert(iraffilename)
Convert the IRAF file name to its Unix equivalent.
Input can be in
directory$file
or$directory/file
format. If'$'
is not found in the input string, it is returned as-is.- Parameters:
iraffilename (str) – Filename in IRAF format.
- Returns:
unixfilename – Filename in Unix format.
- Return type:
- Raises:
AttributeError – Input is not a string.
pysynphot.obsbandpass
This module handle bandpass of observation modes.
- pysynphot.obsbandpass.ObsBandpass(obstring, graphtable=None, comptable=None, component_dict={})
Generate a bandpass object from observation mode.
If the bandpass consists of multiple throughput files (e.g., “acs,hrc,f555w”), then
ObsModeBandpass
is returned. Otherwise, if it consists of a single throughput file (e.g., “johnson,v”), thenTabularSpectralElement
is returned.See Observation Mode and Appendix B: OBSMODE Keywords for more details.
- Parameters:
obstring (str) – Observation mode.
graphtable – See
ObservationMode
.comptable – See
ObservationMode
.component_dict – See
ObservationMode
.
- Returns:
bp
- Return type:
Examples
>>> bp1 = S.ObsBandpass('acs,hrc,f555w') >>> bp2 = S.ObsBandpass('johnson,v')
- class pysynphot.obsbandpass.ObsModeBandpass(ob)
Bases:
CompositeSpectralElement
Bandpass instantiated from an
obsmode
string. Also see Observation Mode, Appendix B: OBSMODE Keywords, and Appendix C: TMG, TMC, and TMT Files.- Parameters:
ob (str) – Observation mode.
- obsmode, name
Same as input
ob
.
- component1, component2
Components and sub-components that belong to the observation mode.
- Type:
CompositeSpectralElement
orSpectralElement
- warnings
To store warnings, which are inherited from all inputs. If they have the same warning keyword, the one from most recently read component is used.
- Type:
- binset
This is set with
bandWave()
, as described in Reference Data.- Type:
array_like
- waveunits
User unit inherited from inputs, where all inputs are required to have the same unit or an exception will be raised.
- Type:
- wave, throughput
Wavelength set in user unit and associated unitless throughput.
- Type:
array_like
- Raises:
NotImplementedError – Inputs have different wavelength units.
TypeError – Both input spectra must be bandpasses.
pysynphot.exceptions.IncompatibleSources – Input spectra have different telescope areas defined.
- pixel_range(waverange, waveunits=None, round='round')
Returns the number of wavelength bins within
waverange
.Note
This calls
pixel_range()
withself.binset
as the first argument.- Parameters:
waverange – See
pixel_range()
.round – See
pixel_range()
.waveunits (str, optional) – The units of the wavelengths given in
waverange
. IfNone
(default), the wavelengths are assumed to be in the units ofself.waveunits
.
- Returns:
num – Number of wavelength bins within
waverange
.- Return type:
- Raises:
pysynphot.exceptions.UndefinedBinset – If
self.binset
isNone
.
- showfiles()
Same as
pysynphot.observationmode.BaseObservationMode.showfiles()
.
- thermback()
Calculate thermal background count rate for
self.obsmode
.Calculation uses
ThermalSpectrum()
to extract thermal component source spectrum inphotlam
per square arcsec. Then this spectrum is integrated and multiplied by detector pixel scale and telescope collecting area to produce a count rate in count/s/pix. This unit is non-standard but used widely by STScI Exposure Time Calculator.Note
Similar to IRAF STSDAS SYNPHOT
thermback
task.- Returns:
ans – Thermal background count rate.
- Return type:
- Raises:
NotImplementedError – Bandpass has no thermal information in graph table.
- wave_range(cenwave, npix, waveunits=None, round='round')
Get the wavelength range covered by the given number of pixels centered on the given wavelength.
Note
This calls
wave_range()
withself.binset
as the first argument.- Parameters:
cenwave – See
wave_range()
.npix – See
wave_range()
.round – See
wave_range()
.waveunits (str, optional) – Wavelength units of
cenwave
and the returned wavelength range. IfNone
(default), the wavelengths are assumed to be in the units ofself.waveunits
.
- Returns:
waverange – The range of wavelengths spanned by
npix
centered oncenwave
.- Return type:
tuple of floats
- Raises:
pysynphot.exceptions.UndefinedBinset – If
self.binset
is None.
- pysynphot.obsbandpass.pixel_range(bins, waverange, round='round')
Returns the number of wavelength bins within
waverange
.- Parameters:
bins (ndarray) – Wavelengths of pixel centers. Must be in the same units as
waverange
.waverange (array_like) – A sequence containing the wavelength range of interest. Only the first and last elements are used. Assumed to be in increasing order. Must be in the same units as
bins
.round ({‘round’, ‘min’, ‘max’,
None
}, optional) –How to deal with pixels at the edges of the wavelength range. All of the options, except
None
, will return an integer number of pixels. Defaults to'round'
.'round'
- Wavelength ends that fall in the middle of a pixel are counted if more than half of the pixel is withinwaverange
. Ends that fall in the center of a pixel are rounded up to the nearest pixel edge.'min'
- Only pixels wholly withinwaverange
are counted.'max'
- End pixels that are withinwaverange
by any margin are counted.None
- The exact number of encompassed pixels, including fractional pixels, is returned.
- Returns:
num – Number of wavelength bins within
waverange
.- Return type:
- Raises:
ValueError – If
round
is not an allowed value.pysynphot.exceptions.OverlapError – If
waverange
exceeds the bounds ofbins
.
- pysynphot.obsbandpass.wave_range(bins, cenwave, npix, round='round')
Get the wavelength range covered by the given number of pixels centered on the given wavelength.
- Parameters:
bins (ndarray) – Wavelengths of pixel centers. Must be in the same units as
cenwave
.cenwave (float) – Central wavelength of range. Must be in the same units as
bins
.npix (int) – Number of pixels in range, centered on
cenwave
.round ({‘round’, ‘min’, ‘max’,
None
}, optional) –How to deal with pixels at the edges of the wavelength range. All of the options, except
None
, will return wavelength ends that correpsonds to pixel edges. Defaults to'round'
.'round'
- A wavelength range is returned such that the ends are pixel edges and the range spans exactlynpix
pixels. Ends that fall in the center of bins are rounded up to the nearest pixel edge.'min'
- The returned wavelength range is shrunk so that it includes an integer number of pixels and the ends fall on pixel edges. May not span exactlynpix
pixels.'max'
- The returned wavelength range is expanded so that it includes an integer number of pixels and the ends fall on pixel edges. May not span exactlynpix
pixels.None
- An exact wavelength range is returned. The wavelength ends returned may not correspond to pixel edges, but will cover exactlynpix
pixels.
- Returns:
waverange – The range of wavelengths spanned by
npix
centered oncenwave
.- Return type:
tuple of floats
- Raises:
ValueError – If
round
is not an allowed value.pysynphot.exceptions.OverlapError – If
cenwave
is not withinbins
, or the returnedwaverange
would exceed the limits ofbins
.
pysynphot.observationmode
This module handles observation modes that are defined in graph tables.
Global Variables
pysynphot.observationmode.rootdir
- Same aspysynphot.locations.rootdir
.pysynphot.observationmode.datadir
- Same aspysynphot.locations.specdir
.pysynphot.observationmode.wavecat
- Same aspysynphot.locations.wavecat
.pysynphot.observationmode.CLEAR
- String to represent a clear filter in an observation mode, i.e., ‘clear’.
- class pysynphot.observationmode.BaseObservationMode(obsmode, method='HSTGraphTable', graphtable=None)
Class that handles the graph table, common to both optical and thermal observation modes. Also see Appendix C: TMG, TMC, and TMT Files.
- Parameters:
- pardict
Stores parameterized keywords and their values. For example,
aper#0.1
would result in{'aper':0.1}
.- Type:
- modes
Individual keywords that make up the observation mode. This includes parameterized ones.
- compnames, thcompnames
Optical and thermal component names based on keyword look-ups. The look-up is done using
GetComponentsFromGT()
.
- components, pixscale
Reserved to be used by sub-classes.
- Type:
- GetFileNames()
Return throughput files of this observation mode.
- Returns:
throughput_filenames
- Return type:
- bandWave()
Return the binned wavelength set most appropriate for the observation mode, as defined by
pysynphot.locations.wavecat
. Also see Reference Data.- Returns:
bandwave
- Return type:
array_like
- showfiles()
Like
GetFileNames()
but print the filenames instead.'clear'
components are not printed.Note
Similar to IRAF STSDAS SYNPHOT
showfiles
task.
- class pysynphot.observationmode.ObservationMode(obsmode, method='HSTGraphTable', graphtable=None, comptable=None, component_dict={})
Bases:
BaseObservationMode
Class to handle optical observation mode.
- Parameters:
obsmode – See
BaseObservationMode
.method – See
BaseObservationMode
.graphtable – See
BaseObservationMode
.comptable (str or
None
) – Component table name. IfNone
, it is taken fromrefs
.component_dict (dict) – Maps component filename to corresponding component object.
- pardict
Stores parameterized keywords and their values. For example,
aper#0.1
would result in{'aper':0.1}
.- Type:
- modes
Individual keywords that make up the observation mode. This includes parameterized ones.
- gtname, ctname
Graph and component table names.
- Type:
- compnames, thcompnames
Optical and thermal component names based on keyword look-ups. The look-up is done using
GetComponentsFromGT()
.
- components
List of component objects. Each object has
throughput_name
(str),throughput
(SpectralElement
), andwaveunits
(Units
) attributes.- Type:
- Raises:
IndexError – Component look-up failed.
- ThermalSpectrum()
Calculate thermal spectrum.
- Returns:
sp – Thermal spectrum in
photlam
.- Return type:
- Raises:
IndexError – Calculation failed.
- Throughput()
Combined throughput from multiplying all the components together.
- Returns:
throughput – Combined throughput.
- Return type:
pysynphot.observation
This module handles an observation and related calculations.
- class pysynphot.observation.Observation(spec, band, binset=None, force=None)
Bases:
CompositeSourceSpectrum
Class to handle an observation. An observation is the end point of a chain of spectral manipulation.
Most
ObsBandpass
objects have a built-inbinset
that is optimized for use with the specified observing mode (also see Wavelength Table). Specifying thebinset
here would override the built-in one.- Parameters:
spec (
SourceSpectrum
) – Source spectrum.band (
SpectralElement
) – Bandpass.binset (array_like or
None
) – Wavelength values to be used for binning when converting to counts. Seeinitbinset()
.force – See
validate_overlap()
.
- spectrum
Same as input
spec
.
- bandpass
Same as input
band
.
- binset
Same as input
binset
.
- component1, component2
Components and sub-components that belong to the observation.
- Type:
- warnings
To store warnings, which are inherited from all inputs. If they have the same warning keyword, the one from most recently read component is used.
- Type:
- isAnalytic
Flag to indicate whether this is an analytic spectrum. This is only
True
if both inputs are analytic.- Type:
- primary_area
Area of the telescope. This is inherited from either of the inputs, if available (not
None
). If inputs have different values, an exception is raised.- Type:
number or
None
- waveunits, fluxunits
User units inherited from source spectrum.
- Type:
- wave, flux
Wavelength set and associated flux in user units. This is the native dataset.
- Type:
array_like
- binwave, binflux
Binned dataset.
- Type:
array_like
- Raises:
pysynphot.exceptions.IncompatibleSources – Input spectra have different telescope areas defined.
- as_spectrum(binned=True)
Reduce the observation to a simple spectrum object.
An observation is a complex object with some restrictions on its capabilities. At times, it would be useful to work with the simulated observation as a simple object that is easier to manipulate and takes up less memory.
- countrate(binned=True, range=None, force=False)
Calculate effective stimulus in count/s. Also see Count Rate and Effective Stimulus.
Note
This is the calculation performed when the ETC invokes
countrate
.- Parameters:
binned (bool) – If
True
(default), use binned data. Otherwise, use native data.range (tuple or
None
) –If not
None
, it must be a sequence with two floating-point elements specifying the wavelength range (inclusive) in the unit ofself.waveunits
in the form of(low, high)
; This is the range over which the integration will be performed. If the specified range does not exactly match a value in the wavelength set:If
binned=True
, the bin containing the range value will be used. This assumesself.binwave
contains bin centers.If
binned=False
, native dataset will be interpolated to the specified values. (Not Implemented.)
force (bool) – If
False
(default), partially overlapping ranges will raise an exception. IfTrue
, a partial overlap will return the calculated value instead. Disjoint ranges raise an exception regardless.
- Returns:
ans – Count rate.
- Return type:
- Raises:
NotImplementedError – Wavelength range is defined for unbinned data.
pysynphot.exceptions.DisjointError – Wavelength range does not overlap with observation.
pysynphot.exceptions.PartialOverlap – Wavelength range only partially overlaps with observation.
- efflam(binned=True)
Calculate effective wavelength of the observation. Calculation is done in the flux unit of
flam
.Note
Similar to IRAF STSDAS SYNPHOT
efflphot
task.
- effstim(fluxunits='photlam')
Compute effective stimulus.
Calculations are done in given flux unit, and wavelengths in Angstrom. Native dataset is used.
- Parameters:
fluxunits (str) – Flux unit.
- Returns:
ans – Effective stimulus.
- Return type:
- Raises:
ValueError – Invalid integrated flux.
- initbinflux()
Calculate binned flux and edges.
Flux is computed by integrating the spectrum on the specified binned wavelength set, using information from the natural wavelength set.
Native wave/flux arrays should be considered samples of a continuous function, but not their binned counterparts. Thus, it makes sense to interpolate
(wave, flux)
but not(binwave, binflux)
.Note
Assumes that the wavelength values in the binned wavelength set are the centers of the bins.
Uses
pysynphot.pysynphot_utils.calcbinflux()
C-extension, if available, for binned flux calculation.
- initbinset(binset=None)
Set
self.binwave
.By default, wavelength values for binning are inherited from bandpass. If the bandpass has no binning information, then source spectrum wavelengths are used. However, if user provides values, then those are used without question.
- Parameters:
binset (array_like or
None
) – Wavelength values to be used for binning when converting to counts.
- pivot(binned=True)
Calculate pivot wavelength of the observation.
Note
This is the calculation performed when ETC invokes
calcphot
.
- pixel_range(waverange, waveunits=None, round='round')
Calculate the number of wavelength bins within given wavelength range.
Note
This calls
pysynphot.obsbandpass.pixel_range()
withself.binwave
as the first argument.- Parameters:
waverange – See
pysynphot.obsbandpass.pixel_range()
.round – See
pysynphot.obsbandpass.pixel_range()
.waveunits (str, optional) – The unit of the wavelength range. If
None
(default), the wavelengths are assumed to be in the units ofself.waveunits
.
- Returns:
num – Number of wavelength bins within
waverange
.- Return type:
- Raises:
pysynphot.exceptions.UndefinedBinset – No binned dataset.
- sample(swave, binned=True, fluxunits='counts')
Sample the observation at the given wavelength. Also see Sampling.
- Parameters:
- Returns:
ans – Sampled flux in given unit.
- Return type:
- Raises:
NotImplementedError – Flux unit is not supported or non-scalar wavelength is given.
ValueError – Given wavelength out of range.
- validate_overlap(force)
Validate that spectrum and bandpass overlap. Warnings are stored in
self.warnings
.- Parameters:
force ({‘extrap’, ‘taper’,
None
}) – IfNone
, it is required that the spectrum and bandpass fully overlap. Partial overlap is allowed if this is set to'extrap'
or'taper'
. Seevalidate_overlap()
.
- wave_range(cenwave, npix, waveunits=None, round='round')
Calculate the wavelength range covered by the given number of pixels, centered on the given wavelength.
Note
This calls
pysynphot.obsbandpass.wave_range()
withself.binwave
as the first argument.- Parameters:
cenwave – See
pysynphot.obsbandpass.wave_range()
.npix – See
pysynphot.obsbandpass.wave_range()
.round – See
pysynphot.obsbandpass.wave_range()
.waveunits (str, optional) – Wavelength unit of the given and the returned wavelength values. If
None
(default), the wavelengths are assumed to be in the unit ofself.waveunits
.
- Returns:
waverange – The range of wavelengths spanned by
npix
centered oncenwave
.- Return type:
tuple of floats
- Raises:
pysynphot.exceptions.UndefinedBinset – No binned dataset.
- writefits(fname, clobber=True, trimzero=True, binned=True, hkeys=None)
Like
pysynphot.spectrum.SourceSpectrum.writefits()
but withbinned=True
as default.
- pysynphot.observation.check_overlap(a, b)
Check for wavelength overlap between two spectra.
Note
Generalized from
pysynphot.spectrum.SpectralElement.check_overlap()
.- Parameters:
a (
SourceSpectrum
orSpectralElement
) – Typically a source spectrum, spectral element, observation, or bandpass from observation mode.b (
SourceSpectrum
orSpectralElement
) – Typically a source spectrum, spectral element, observation, or bandpass from observation mode.
- Returns:
result – Full, partial, or no overlap.
- Return type:
{‘full’, ‘partial’, ‘none’}
- Raises:
AttributeError – Given spectrum does not have flux or throughput.
- pysynphot.observation.validate_overlap(comp1, comp2, force)
Validate the overlap between the wavelength sets of the two given components.
- Parameters:
comp1 (
SourceSpectrum
orSpectralElement
) – Source spectrum and bandpass of an observation.comp2 (
SourceSpectrum
orSpectralElement
) – Source spectrum and bandpass of an observation.force ({‘extrap’, ‘taper’,
None
}) – If notNone
, the components may be adjusted by extrapolation or tapering.
- Returns:
comp1, comp2 – Same as inputs. However,
comp1
might be tapered if that option is selected.warnings (dict) – Maps warning keyword to its description.
- Raises:
KeyError – Invalid
force
.pysynphot.exceptions.DisjointError – No overlap detected when
force
isNone
.pysynphot.exceptions.PartialOverlap – Partial overlap detected when
force
isNone
.
pysynphot.planck
This module handles Planck’s law of blackbody radiation.
Global Variables
pysynphot.planck.H
- Planck’s constant in CGS units.pysynphot.planck.HS
- Planck’s constant in SI units.pysynphot.planck.C
- Speed of light in SI units.pysynphot.planck.K
- Boltzmann constant in SI units.
These are used in calculations to prevent floating point overflow, as defined
in IRAF STSDAS SYNPHOT bbfunc
task:
pysynphot.planck.LOWER
pysynphot.planck.UPPER
These are constants used in llam_SI()
:
pysynphot.planck.C1
- Power \(\times\) unit area per steradian.pysynphot.planck.C2
This is used in bb_photlam_arcsec()
:
pysynphot.planck.F
- Factor for conversion from \(\mathrm{m}^{2} \; \mathrm{sr}^{-1} \; \mathrm{m}^{-1}\) to \(\mathrm{cm}^{2} \; \mathrm{arcsec}^{-2} \; \AA^{-1}\).
- pysynphot.planck.bb_photlam_arcsec(wave, temperature)
Evaluate Planck’s law in
photlam
per square arcsec.Note
Uses
llam_SI()
for calculation, and then converts SI units back to CGS.- Parameters:
wave (array_like) – Wavelength values in Angstrom.
temperature (float) – Blackbody temperature in Kelvin.
- Returns:
result – Blackbody radiation in
photlam
per square arcsec.- Return type:
array_like
- pysynphot.planck.bbfunc(wave, temperature)
Evaluate Planck’s law in
photlam
(per steradian).Note
Adapted from IRAF STSDAS SYNPHOT
bbfunc
task.- Parameters:
wave (array_like) – Wavelength values in Angstrom.
temperature (float) – Blackbody temperature in Kelvin.
- Returns:
result – Blackbody radiation in
photlam
per steradian.- Return type:
array_like
- pysynphot.planck.llam_SI(wave, temperature)
Like
bbfunc()
but in SI units.Note
Adapted from SSP code by Dr. Anand Sivaramakrishnan in IRAF STSDAS SYNPHOT.
- Parameters:
wave (array_like) – Wavelength values in meters.
temperature (float) – Blackbody temperature in Kelvin.
- Returns:
result – Blackbody radiation in SI units.
- Return type:
array_like
pysynphot.reddening
This module handles reddening laws and extinction calculations.
- pysynphot.reddening.Extinction(extval, name=None)
Generate extinction curve to be used with spectra.
By default,
reddening()
is used to generate the extinction curve. If a deprecated reddening law is given, thenDeprecatedExtinction
is used instead.Note
Reddening laws are cached in
pysynphot.Cache.RedLaws
for better performance. Repeated calls to the same reddening law here returns the cached result.- Parameters:
extval (float) – Value of \(E(B-V)\) in magnitudes.
name (str or
None
) – Name of reddening law (seeprint_red_laws()
). IfNone
(default), the average Milky Way extinction ('mwavg'
) will be used.
- Returns:
ext – Extinction curve.
- Return type:
- Raises:
ValueError – Invalid reddening law.
Examples
>>> ext = S.Extinction(0.3, 'mwavg')
- pysynphot.reddening.print_red_laws()
Print available extinction laws to screen.
Available extinction laws are extracted from
pysynphot.locations.EXTDIR
. The printed names may be used withExtinction()
to retrieve available reddening laws.Examples
>>> S.reddening.print_red_laws() name reference -------- -------------------------------------------------------------- None Cardelli, Clayton, & Mathis (1989, ApJ, 345, 245) R_V = 3.10. gal3 Cardelli, Clayton, & Mathis (1989, ApJ, 345, 245) R_V = 3.10. lmc30dor Gordon et al. (2003, ApJ, 594, 279) R_V = 2.76. lmcavg Gordon et al. (2003, ApJ, 594, 279) R_V = 3.41. mwavg Cardelli, Clayton, & Mathis (1989, ApJ, 345, 245) R_V = 3.10. mwdense Cardelli, Clayton, & Mathis (1989, ApJ, 345, 245) R_V = 5.00. mwrv21 Cardelli, Clayton, & Mathis (1989, ApJ, 345, 245) R_V = 2.1. mwrv4 Cardelli, Clayton, & Mathis (1989, ApJ, 345, 245) R_V = 4.0. smcbar Gordon et al. (2003, ApJ, 594, 279) R_V=2.74. xgalsb Calzetti et al. (2000. ApJ, 533, 682)
- class pysynphot.reddening.CustomRedLaw(wave=None, waveunits='InverseMicrons', Avscaled=None, name='Unknown Reddening Law', litref=None)
Class to handle reddening law.
- Parameters:
- wave, waveunits, name, litref
Same as inputs.
- obscuration
Same as input
Avscaled
.
- reddening(extval)
Compute the reddening for the given extinction.
\[ \begin{align}\begin{aligned}A(V) = R(V) \; \times \; E(B-V)\\\mathrm{THRU} = 10^{-0.4 \; A(V)}\end{aligned}\end{align} \]Note
self.litref
is passed intoans.citation
.- Parameters:
extval (float) – Value of \(E(B-V)\) in magnitudes.
- Returns:
ans – Extinction curve to apply to a source spectrum.
- Return type:
- class pysynphot.reddening.RedLaw(filename)
Bases:
CustomRedLaw
CustomRedLaw
from a FITS file.Table must be in EXT 1 and contains the following columns:
WAVELENGTH
Av/E(B-V)
Wavelength unit is extracted from
TUNIT1
keyword in EXT 1 header. The primary header (EXT 0) must haveSHORTNM
andLITREF
keywords.- Parameters:
filename (str) – FITS table filename.
- wave
Wavelength values from the
WAVELENGTH
column in EXT 1.- Type:
array_like
- obscuration
Values from the
Av/E(B-V)
column in EXT 1.
pysynphot.refs
This module handles constants and look-up tables used in calculations.
Global Variables
pysynphot.refs._default_waveset
- Default wavelength set to use if no instrument-specific values found.pysynphot.refs._default_waveset_str
- Description of the default wavelength set above.pysynphot.refs.PRIMARY_AREA
- Telescope collecting area, i.e., the primary mirror, in \(\mathrm{cm}^{2}\). The value for HST is 45238.93416.
These are used in observationmode
to look up throughput files for
a given bandpass:
pysynphot.refs.GRAPHTABLE
pysynphot.refs.GRAPHDICT
pysynphot.refs.COMPTABLE
pysynphot.refs.COMPDICT
pysynphot.refs.THERMTABLE
pysynphot.refs.THERMDICT
- pysynphot.refs.getref()
Current default values for graph and component tables, primary area, and wavelength set.
Note
Also see
setref()
.- Returns:
ans – Mapping of parameter names to their current values.
- Return type:
- pysynphot.refs.set_default_waveset(minwave=500, maxwave=26000, num=10000, delta=None, log=True)
Set the default wavelength set,
pysynphot.refs._default_waveset
.- Parameters:
minwave (float, optional) – The start (inclusive) and end (exclusive) points of the wavelength set. Values should be given in linear space regardless of
log
.maxwave (float, optional) – The start (inclusive) and end (exclusive) points of the wavelength set. Values should be given in linear space regardless of
log
.num (int, optional) – The number of elements in the wavelength set. Only used if
delta=None
.delta (float, optional) – Delta between values in the wavelength set. If
log=True
, this defines wavelegth spacing in log space.log (bool, optional) – Determines whether the wavelength set is evenly spaced in log or linear space.
- pysynphot.refs.setref(graphtable=None, comptable=None, thermtable=None, area=None, waveset=None)
Set default graph and component tables, primary area, and wavelength set.
This is similar to setting
refdata
in IRAF STSDAS SYNPHOT. If all parameters set toNone
, they are reverted to software default. If any of the parameters are notNone
, they are set to desired values while the rest (if any) remain at current setting.- Parameters:
graphtable (str or
None
) – Graph, component, and thermal table names, respectively, forobservationmode
throughput look-up. Do not use “*” wildcard.comptable (str or
None
) – Graph, component, and thermal table names, respectively, forobservationmode
throughput look-up. Do not use “*” wildcard.thermtable (str or
None
) – Graph, component, and thermal table names, respectively, forobservationmode
throughput look-up. Do not use “*” wildcard.area (float or
None
) – Telescope collecting area, i.e., the primary mirror, in \(\mathrm{cm}^{2}\).waveset (tuple or
None
) –- Parameters for
set_default_waveset()
as follow: (minwave, maxwave, num)
- This assumes log scale.(minwave, maxwave, num, 'log')
(minwave, maxwave, num, 'linear')
- Parameters for
- Raises:
ValueError – Invalid
waveset
parameters.
pysynphot.renorm
This module handles normalization of source spectrum flux.
- pysynphot.renorm.DefineStdSpectraForUnits()
Define
StdSpectrum
attribute for all the supported Flux Units.This is automatically done on module import. The attribute stores the source spectrum necessary for normalization in the corresponding flux unit.
For
photlam
,photnu
,flam
,fnu
, Jy, and mJy, the spectrum is flat in the respective units with flux value of 1.For counts and
obmag
, it is flat in the unit of counts with flux value of \(1 / N\), where \(N\) is the size of default wavelength set (seerefs
).For
abmag
andstmag
, it is flat in the respective units with flux value of 0 mag. That is equivalent to \(3.63 \times 10^{-20}\)fnu
and \(3.63 \times 10^{-9}\)flam
, respectively.For
vegamag
, it is simply Vega.
- pysynphot.renorm.StdRenorm(spectrum, band, RNval, RNunitstring, force=False)
This is used by
SourceSpectrum
for renormalization.- Parameters:
spectrum (
SourceSpectrum
) – Spectrum to renormalize.band – See
renorm()
.RNval – See
renorm()
.RNunitstring – See
renorm()
.force – See
renorm()
.
- Returns:
newsp – Renormalized spectrum.
- Return type:
pysynphot.spark
SPARK 0.6.1 is obtained from John Aycock (1998-2000). It is the underlying engine used by pysynphot.spparser.
pysynphot.spectrum
This module contains the basis for all spectra classes, including source spectra and bandpasses.
It also pre-loads the built-in Vega spectrum to
pysynphot.spectrum.Vega
.
- class pysynphot.spectrum.Integrator
Integrator engine, which is the base class for
SourceSpectrum
andSpectralElement
.- trapezoidIntegration(x, y)
Perform trapezoid integration.
- Parameters:
x (array_like) – Wavelength set.
y (array_like) – Integrand. For example, throughput or throughput multiplied by wavelength.
- Returns:
sum – Integrated sum.
- Return type:
- validate_fluxtable()
Check for non-negative fluxes. If found, the negative flux values are set to zero, and a warning is printed to screen. This check is not done if flux unit is a magnitude because negative magnitude values are legal.
- validate_wavetable()
Enforce monotonic, ascending wavelength array with no zero or negative values.
- Raises:
pysynphot.exceptions.DuplicateWavelength – Wavelength array contains duplicate entries.
pysynphot.exceptions.UnsortedWavelength – Wavelength array is not monotonic ascending or descending.
pysynphot.exceptions.ZeroWavelength – Wavelength array has zero or negative value(s).
- class pysynphot.spectrum.SourceSpectrum
Bases:
Integrator
This is the base class for all source spectra.
- addmag(magval)
Add a scalar magnitude to existing flux values.
\[\mathrm{flux}_{\mathrm{new}} = 10^{-0.4 \; \mathrm{magval}} \; \mathrm{flux}\]- Parameters:
magval (number) – Magnitude value.
- Returns:
sp – New source spectrum with adjusted flux values.
- Return type:
- Raises:
TypeError – Magnitude value is not a scalar number.
- convert(targetunits)
Set new user unit, for either wavelength or flux.
This effectively converts the spectrum wavelength or flux to given unit. Note that actual data are always kept in internal units (Angstrom and
photlam
), and only converted to user units bygetArrays()
during actual computation. User units are stored inself.waveunits
andself.fluxunits
.
- property flux
Flux property.
- getArrays()
Return wavelength and flux arrays in user units.
- Returns:
wave (array_like) – Wavelength array in
self.waveunits
.flux (array_like) – Flux array in
self.fluxunits
. When necessary,self.primary_area
is used for unit conversion.
- integrate(fluxunits='photlam')
Integrate the flux in given unit.
Integration is done using
trapezoidIntegration()
withx=wave
andy=flux
, where flux has been convert to given unit first.\[\mathrm{result} = \int F_{\lambda} d\lambda\]
- redshift(z)
Apply redshift to the spectrum.
Redshifted spectrum is never analytic even if the input spectrum is. Output units are always Angstrom and PHOTLAM regardless of user units.
- Parameters:
z (number) – Redshift value.
- Returns:
copy – Redshifted spectrum.
- Return type:
- renorm(RNval, RNUnits, band, force=False)
Renormalize the spectrum to the specified value, unit, and bandpass.
This wraps
StdRenorm()
for convenience. Basically, the spectrum is multiplied by a numeric factor so that the total integrated flux will be the given value in the given unit in the given bandpass.When
force=False
, if spectrum is not fully defined within the given bandpass, but the overlap is at least 99%, a warning is printed to screen andself.warnings['PartialRenorm']
is set toTrue
.- Parameters:
RNval (number) – Flux value for renormalization.
band (
SpectralElement
) – Bandpass thatRNval
is based on.force (bool) – Force renormalization regardless of overlap status with given bandpass. If
True
, overlap check is skipped. Default isFalse
.
- Returns:
newsp – Renormalized spectrum.
- Return type:
- Raises:
ValueError – Integrated flux is zero, negative, NaN, or infinite.
pysynphot.exceptions.DisjointError – Spectrum and bandpass are disjoint.
pysynphot.exceptions.OverlapError – Spectrum and bandpass do not fully overlap.
- sample(wave, interp=True)
Sample the spectrum at given wavelength(s).
This method has two behaviors:
When
interp=True
, wavelength(s) must be provided as a Numpy array. Interpolation is done in internal units (Angstrom andphotlam
).When
interp=False
, wavelength must be a scalar number. The flux that corresponds to the closest matching wavelength value is returned. This option should only be used for sampling binned data inObservation
.
- Parameters:
- Returns:
ans – Sampled flux in user unit.
- Return type:
number or array_like
- Raises:
NotImplementedError – Non-scalar wavelength set provided when interpolation is not allowed.
- validate_units()
Ensure that wavelenth and flux units belong to the correct classes.
- property wave
Wavelength property.
- writefits(filename, clobber=True, trimzero=True, binned=False, precision=None, hkeys=None)
Write the spectrum to a FITS table.
Primary header in EXT 0.
FILENAME
,ORIGIN
, and any extra keyword(s) fromhkeys
will also be added.Table header and data are in EXT 1. The table has 2 columns, i.e.,
WAVELENGTH
andFLUX
. Data are stored in user units. Its header also will have these additional keywords:EXPR
- Description of the spectrum.TDISP1
andTDISP2
- Columns display format, always “G15.7”.GRFTABLE
andCMPTABLE
- Graph and component table names to use with associated observation mode. These are only added if applicable.
If data is already double-precision but user explicitly set output precision to single,
pysynphot.spectrum.syn_epsilon
defines the allowed minimum wavelength separation. This limit (\(3.2 \times 10^{-4}\)) was taken from IRAF STSDAS SYNPHOT FAQ. Values equal or smaller than this limit are considered as the same, and duplicates are ignored, resulting in data loss. In the way that this comparison is coded, when such precision clash happens, even when no duplicates are detected, the last row is always omitted (also data loss). Therefore, it is not recommended for user to force single-precision when the data is in double-precision.- Parameters:
filename (str) – Output filename.
trimzero (bool) – Trim off duplicate rows with flux values of zero from both ends of the spectrum. This keeps one row of zero-flux at each end, if it exists; However, it does not add a zero-flux row if it does not. Default is
True
.binned (bool) – Write
self.binwave
andself.binflux
(binned) dataset, instead ofself.wave
andself.flux
(native). Using this option when object does not have binned data will cause an exception to be raised. Default isFalse
.precision ({‘s’, ‘d’,
None
}) – Write data out in single ('s'
) or double ('d'
) precision. Default isNone
, which will enforce native precision fromself.flux
.hkeys (dict) – Additional keyword(s) to be added to primary FITS header, in the format of
{keyword:(value,comment)}
.
- class pysynphot.spectrum.TabularSourceSpectrum(filename=None, fluxname=None, keepneg=False)
Bases:
SourceSpectrum
Base class for
ArraySourceSpectrum
andFileSourceSpectrum
.- Parameters:
filename (str or
None
) – File with spectral data (can be ASCII or FITS). If notNone
, data will be loaded from file at initialization.fluxname (str or
None
) – Column name containing flux data. This is only used if filename is given and is of FITS format.keepneg (bool) – Keep negative flux values instead of setting them to zero with a warning. Default is
False
.
- filename, name
Same as input.
- waveunits, fluxunits
User units for wavelength and flux.
- Type:
- wave, flux
Wavelength set and associated flux in user units.
- Type:
array_like
- GetWaveSet()
Return the wavelength set for the spectrum.
- Returns:
waveset – Wavelength set (a copy of the internal wavelength table).
- Return type:
array_like
- resample(resampledWaveTab)
Resample the spectrum for the given wavelength set.
Given wavelength array must be monotonically increasing or decreasing. Flux interpolation is done using
numpy.interp()
.- Parameters:
resampledWaveTab (array_like) – Wavelength set for resampling.
- Returns:
resampled – Resampled spectrum.
- Return type:
- taper()
Taper the spectrum by adding zero flux to each end. This is similar to
SpectralElement.taper()
.There is no check to see if the spectrum is already tapered. Hence, calling this on a tapered spectrum will result in multiple zero-flux entries at both ends.
The wavelengths to use for the new first and last points are calculated by using the same ratio as for the two interior points used at each end.
- Returns:
OutSpec – Tapered spectrum.
- Return type:
- class pysynphot.spectrum.ArraySourceSpectrum(wave=None, flux=None, waveunits='angstrom', fluxunits='photlam', name='UnnamedArraySpectrum', keepneg=False)
Bases:
TabularSourceSpectrum
Class to handle source spectrum from arrays.
- Parameters:
wave (array_like) – Wavelength and flux arrays.
flux (array_like) – Wavelength and flux arrays.
waveunits (str) – Wavelength and flux units, as accepted by
Units
. Defaults are Angstrom andphotlam
.fluxunits (str) – Wavelength and flux units, as accepted by
Units
. Defaults are Angstrom andphotlam
.name (str) – Description of the spectrum. Default is “UnnamedArraySpectrum”.
keepneg (bool) – Keep negative flux values instead of setting them to zero with a warning. Default is
False
.
- name
Same as input.
- waveunits, fluxunits
User units for wavelength and flux.
- Type:
- wave, flux
Wavelength set and associated flux in user units.
- Type:
array_like
- Raises:
ValueError – Mismatched wavelength and flux arrays.
- class pysynphot.spectrum.FileSourceSpectrum(filename, fluxname=None, keepneg=False)
Bases:
TabularSourceSpectrum
Class to handle source spectrum loaded from ASCII or FITS table. Also see File I/O.
- Parameters:
- fheader
For FITS file, this contains headers from both extensions 0 and 1. If the extensions have the same keyword, the one from the latter is used.
- Type:
- waveunits, fluxunits
User units for wavelength and flux.
- Type:
- wave, flux
Wavelength set and associated flux in user units.
- Type:
array_like
- class pysynphot.spectrum.AnalyticSpectrum(waveunits='angstrom', fluxunits='photlam')
Bases:
SourceSpectrum
Base class for analytic source spectrum. This includes
BlackBody
,FlatSpectrum
,GaussianSource
, andPowerlaw
.- Parameters:
- waveunits, fluxunits
User units for wavelength and flux.
- Type:
- wave, flux
Wavelength set and associated flux in user units.
- Type:
array_like
- GetWaveSet()
Return the wavelength set for the spectrum.
- Returns:
waveset – Wavelength set (a copy of the default wavelength table).
- Return type:
array_like
- class pysynphot.spectrum.GaussianSource(flux, center, fwhm, waveunits='angstrom', fluxunits='flam')
Bases:
AnalyticSpectrum
Class to handle a Gaussian source.
- Parameters:
flux (float) – Total flux under the Gaussian curve, in given flux unit.
center (float) – Central wavelength of the Gaussian curve, in given wavelength unit.
fwhm (float) – FWHM of the Gaussian curve, in given wavelength unit.
waveunits (str) – Wavelength and flux units, as accepted by
Units
. Defaults are Angstrom andflam
.fluxunits (str) – Wavelength and flux units, as accepted by
Units
. Defaults are Angstrom andflam
.
- total_flux
Same as input
flux
.
- center, fwhm
Same as inputs.
- sigma, factor
These are \(\sigma\) and \(A\) as defined in Gaussian Emission.
- Type:
- waveunits, fluxunits
User units for wavelength and flux.
- Type:
- wave, flux
Wavelength set and associated flux in user units.
- Type:
array_like
- GetWaveSet()
Return the wavelength set that optimally samples the Gaussian curve. It has 101 values, as defined below:
\[ \begin{align}\begin{aligned}x_{\mathrm{first,last}} = x_{0} \; \pm \; 5 \; \sigma\\\delta x = 0.1 \; \sigma\end{aligned}\end{align} \]- Returns:
waveset – Wavelength set in internal unit.
- Return type:
array_like
- class pysynphot.spectrum.FlatSpectrum(fluxdensity, waveunits='angstrom', fluxunits='photlam')
Bases:
AnalyticSpectrum
Class to handle a flat source spectrum.
- Parameters:
- waveunits, fluxunits
User units for wavelength and flux.
- Type:
- wave, flux
Wavelength set and associated flux in user units.
- Type:
array_like
- redshift(z)
Apply redshift to the flat spectrum.
Unlike
SourceSpectrum.redshift()
, the redshifted spectrum remains an analytic flat source.- Parameters:
z (number) – Redshift value.
- Returns:
ans
- Return type:
- class pysynphot.spectrum.Powerlaw(refwave, index, waveunits='angstrom', fluxunits='photlam')
Bases:
AnalyticSpectrum
Class to handle a power-law source spectrum.
- Parameters:
- waveunits, fluxunits
User units for wavelength and flux.
- Type:
- wave, flux
Wavelength set and associated flux in user units.
- Type:
array_like
- class pysynphot.spectrum.BlackBody(temperature)
Bases:
AnalyticSpectrum
Class to handle a blackbody source.
Flux is evaluated with
bbfunc()
and normalized withpysynphot.spectrum.RENORM
, which is:\[\mathrm{RENORM} = \pi \; (\frac{R_{\odot}}{1 \; \mathrm{kpc}})^{2}\]- Parameters:
temperature (number) – Blackbody temperature in Kelvin.
- temperature
Same as input.
- waveunits, fluxunits
User units for wavelength and flux.
- Type:
- wave, flux
Wavelength set and associated flux in user units.
- Type:
array_like
- class pysynphot.spectrum.CompositeSourceSpectrum(source1, source2, operation)
Bases:
SourceSpectrum
Class to handle composite spectrum involving source spectra.
- Parameters:
source1 (
SourceSpectrum
orSpectralElement
) – One or both of the inputs must be source spectrum.source2 (
SourceSpectrum
orSpectralElement
) – One or both of the inputs must be source spectrum.operation ({'add', 'multiply'}) – Math operation to perform.
- component1, component2
Same as input
source1
andsource2
.
- operation
Same as input.
- warnings
To store warnings, which are inherited from both input sources. If inputs have the same warning keyword, the one from
source2
is used.- Type:
- isAnalytic
Flag to indicate whether this is an analytic spectrum. This is only
True
if both inputs are analytic.- Type:
- primary_area
Area of the telescope. This is inherited from either of the inputs, if available (not
None
). If inputs have different values, an exception is raised.- Type:
number or
None
- waveunits, fluxunits
User units inherited from
source1
(if available) orsource2
(if not).- Type:
- wave, flux
Wavelength set and associated flux in user units.
- Type:
array_like
- Raises:
pysynphot.exceptions.IncompatibleSources – Input spectra have different telescope areas defined.
- GetWaveSet()
Obtain the wavelength set for the composite spectrum. This is done by using
MergeWaveSets()
to form a union of wavelength sets from its components.- Returns:
waveset – Composite wavelength set.
- Return type:
array_like
- complist()
Return a list of all components and sub-components. This is for use with
__iter__()
.
- tabulate()
Return a simplified version of the spectrum.
Composite spectrum can be overly complicated when it has too many components and sub-components. This method copies the following into a simple (tabulated) source spectrum:
Name
Wavelength array and unit
Flux array and unit
- Returns:
sp – Tabulated source spectrum.
- Return type:
- class pysynphot.spectrum.SpectralElement
Bases:
Integrator
This is the base class for all bandpasses and spectral elements (e.g., filter and detector response curves).
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- GetThroughput()
Obtain throughput for the spectrum.
- Returns:
throughput – Throughput values.
- Return type:
array_like
- GetWaveSet()
Obtain the wavelength set for the spectrum.
- Returns:
wave – Wavelength set in internal unit.
- Return type:
array_like
- avgwave()
Calculate Bandpass Average Wavelength.
- Returns:
ans – Average wavelength.
- Return type:
- check_overlap(other)
Check overlap with another spectrum. Also see Overlap Checks.
This checks whether the wavelength set of the given spectrum is defined everywhere within
self
. Wavelength values where throughput is zero are excluded from the check. Typical use case is for checking whether a source spectrum is fully defined over the range of a bandpass. This check is asymmetric in the sense that ifself
is fully defined within the given spectrum, but not the other way around, it will still only return “partial”. If the given spectrum is analytic, the result is always “full”.Example of full overlap:
|---------- other ----------| |------ self ------|
Examples of partial overlap:
|---------- self ----------| |------ other ------| |---- other ----| |---- self ----| |---- self ----| |---- other ----|
Examples of no overlap:
|---- self ----| |---- other ----| |---- other ----| |---- self ----|
- Parameters:
other (
SourceSpectrum
orSpectralElement
) – The other spectrum.- Returns:
ans – Overlap status.
- Return type:
{‘full’, ‘partial’, ‘none’}
- check_sig(other)
Check overlap insignificance with another spectrum. Also see Overlap Checks.
Note
Only use when
check_overlap()
returns “partial”.- Parameters:
other (
SourceSpectrum
orSpectralElement
) – The other spectrum.- Returns:
ans –
True
means the lack of overlap is insignificant (i.e., okay to proceed).- Return type:
- convert(targetunits)
Set new user unit, for wavelength only.
This effectively converts the spectrum wavelength to given unit. Note that actual data are always kept in internal unit (Angstrom), and only converted to user unit by
GetWaveSet()
during actual computation. User unit is stored inself.waveunits
. Throughput is unitless and cannot be converted.
- efficiency()
Calculate Bandpass Dimensionless Efficiency.
- Returns:
ans – Bandpass dimensionless efficiency.
- Return type:
- equivwidth()
Calculate Bandpass Equivalent Width. This basically just calls
integrate()
.- Returns:
ans – Bandpass equivalent width.
- Return type:
- integrate(wave=None)
Integrate the throughput over the specified wavelength set. If no wavelength set is specified, the built-in one is used.
Integration is done using
trapezoidIntegration()
withx=wave
andy=throughput
. Also see Bandpass Equivalent Width.
- photbw(floor=0)
Calculate Bandpass RMS Band Width (SYNPHOT).
Note
For backward-compatibility with IRAF STSDAS SYNPHOT only.
- Parameters:
floor (float) – Same as
rmswidth()
.- Returns:
ans – RMS band width (deprecated).
- Return type:
- pivot(binned=False)
Calculate Pivot Wavelength.
- Parameters:
binned (bool) – This is reserved for use by
Observation
. IfTrue
, binned wavelength set is used. Default isFalse
.- Returns:
ans – Pivot wavelength.
- Return type:
- Raises:
AttributeError – Binned wavelength set requested but not found.
- rectwidth()
Calculate Bandpass Rectangular Width.
- Returns:
ans – Bandpass rectangular width.
- Return type:
- resample(resampledWaveTab)
Resample the spectrum for the given wavelength set.
Given wavelength array must be monotonically increasing or decreasing. Throughput interpolation is done using
numpy.interp()
.- Parameters:
resampledWaveTab (array_like) – Wavelength set for resampling.
- Returns:
resampled – Resampled spectrum.
- Return type:
- rmswidth(floor=0)
Calculate Bandpass RMS Band Width (Koornneef).
- sample(wave)
Sample the spectrum.
This uses
resample()
to do the actual computation.- Parameters:
wave (number or array_like) – Wavelength set for sampling, given in user unit.
- Returns:
throughput – Sampled throughput.
- Return type:
number or array_like
- taper()
Taper the spectrum by adding zero throughput to each end. This is similar to
TabularSourceSpectrum.taper()
.There is no check to see if the spectrum is already tapered. Hence, calling this on a tapered spectrum will result in multiple zero-throughput entries at both ends.
The wavelengths to use for the new first and last points are calculated by using the same ratio as for the two interior points used at each end.
- Returns:
OutElement – Tapered spectrum.
- Return type:
- property throughput
Throughput property.
- unit_response()
Calculate Bandpass Unit Response.
Warning
Result is correct only if
self.waveunits
is in Angstrom.- Returns:
ans – Bandpass unit response.
- Return type:
- validate_units()
Ensure that wavelenth unit belongs to the correct class. There is no check for throughput because it is unitless.
- property wave
Wavelength property.
- writefits(filename, clobber=True, trimzero=True, precision=None, hkeys=None)
Write the spectrum to a FITS table.
Primary header in EXT 0.
FILENAME
,ORIGIN
, and any extra keyword(s) fromhkeys
will also be added.Table header and data are in EXT 1. The table has 2 columns, i.e.,
WAVELENGTH
andTHROUGHPUT
. Wavelength data are stored in user unit. Its header also will have these additional keywords:EXPR
- Description of the spectrum.TDISP1
andTDISP2
- Columns display format, always “G15.7”.GRFTABLE
andCMPTABLE
- Graph and component table names to use with associated observation mode. These are only added if applicable.
If data is already double-precision but user explicitly set output precision to single,
pysynphot.spectrum.syn_epsilon
defines the allowed minimum wavelength separation. This limit (\(3.2 \times 10^{-4}\)) was taken from IRAF STSDAS SYNPHOT FAQ. Values equal or smaller than this limit are considered as the same, and duplicates are ignored, resulting in data loss. In the way that this comparison is coded, when such precision clash happens, even when no duplicates are detected, the last row is always omitted (also data loss). Therefore, it is not recommended for user to force single-precision when the data is in double-precision.- Parameters:
filename (str) – Output filename.
trimzero (bool) – Trim off duplicate rows with flux values of zero from both ends of the spectrum. This keeps one row of zero-flux at each end, if it exists; However, it does not add a zero-flux row if it does not. Default is
True
.precision ({‘s’, ‘d’,
None
}) – Write data out in single ('s'
) or double ('d'
) precision. Default isNone
, which will enforce native precision fromself.throughput
.hkeys (dict) – Additional keyword(s) to be added to primary FITS header, in the format of
{keyword:(value,comment)}
.
- class pysynphot.spectrum.TabularSpectralElement(fileName=None, thrucol='throughput')
Bases:
SpectralElement
Base class for
ArraySpectralElement
andFileSpectralElement
.- Parameters:
- name
Same as input
fileName
.
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- wave, throughput
Wavelength set in user unit and associated unitless throughput.
- Type:
array_like
- class pysynphot.spectrum.ArraySpectralElement(wave=None, throughput=None, waveunits='angstrom', name='UnnamedArrayBandpass')
Bases:
TabularSpectralElement
Class to handle bandpass from arrays.
- Parameters:
- name
Same as input.
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- wave, throughput
Wavelength set in user unit and associated unitless throughput.
- Type:
array_like
- Raises:
ValueError – Mismatched wavelength and throughput arrays.
- class pysynphot.spectrum.FileSpectralElement(filename, thrucol=None)
Bases:
TabularSpectralElement
Class to handle bandpass loaded from ASCII or FITS table. Also see File I/O.
- Parameters:
- fheader
For FITS file, this contains headers from both extensions 0 and 1. If the extensions have the same keyword, the one from the latter is used.
- Type:
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- wave, throughput
Wavelength set in user unit and associated unitless throughput.
- Type:
array_like
- class pysynphot.spectrum.InterpolatedSpectralElement(fileName, wavelength)
Bases:
SpectralElement
Class to handle parameterized keyword in an observation mode.
- Parameters:
fileName (str) – Filename followed by a column name specification between square brackets. For example: “mythru_syn.fits[fr388n#]”
wavelength (number) – Desired value to interpolate to. This is not restricted to wavelength, but rather whatever parameter the file is parameterized for.
- interpval
Same as input
wavelength
.
- warnings
To store warnings. When extrapolation is not allowed but a default throughput column is present and used,
warnings['DefaultThroughput']
is set toTrue
.- Type:
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- throughputunits
This is only to inform user that throughput is unitless.
- Type:
‘none’
- wave, throughput
Wavelength set in user unit and associated unitless throughput.
- Type:
array_like
- Raises:
Exception – File does not have columns needed for interpolation.
pysynphot.exceptions.ExtrapolationNotAllowed – Extrapolation is not allowed and no default throughput column found.
- class pysynphot.spectrum.UniformTransmission(value, waveunits='angstrom')
Bases:
SpectralElement
Class to handle a uniform bandpass.
- Parameters:
- value
Same as input.
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- wave, throughput
Wavelength set in user unit and associated unitless throughput.
- Type:
array_like
- GetWaveSet()
Obtain wavelength set for the spectrum.
- Returns:
waveset – Due to the nature of uniform transmission, this is always undefined.
- Return type:
- property wave
waveset
for uniform transmission.
- writefits(*args, **kwargs)
Write to file using default waveset.
- class pysynphot.spectrum.Box(center, width, waveunits=None)
Bases:
SpectralElement
Class to handle a box-shaped bandpass.
- Parameters:
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- wave, throughput
Wavelength set in user unit and associated unitless throughput.
- Type:
array_like
- resample(resampledWaveTab)
Resample the spectrum for the given wavelength set.
Given wavelength array must be monotonically increasing or decreasing.
- Parameters:
resampledWaveTab (array_like) – Wavelength set for resampling.
- Returns:
resampled – Resampled spectrum. This is no longer a real
Box
spectrum.- Return type:
- sample(wavelength)
Input wavelengths assumed to be in user unit.
- class pysynphot.spectrum.ThermalSpectralElement(fileName)
Bases:
TabularSpectralElement
Class to handle spectral element with thermal properties read from a FITS table.
Note
This class does not know how to apply itself to an existing beam. Its emissivity is handled by
ThermalSpectrum()
.- Parameters:
fileName (str) – Filename of the thermal emissivity table.
- name
Same as input
fileName
.
- temperature
Default temperature in Kelvin from header.
- Type:
number
- beamFillFactor
Beam filling factor from header.
- Type:
number
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- throughputunits
This is only to inform user that throughput is unitless.
- Type:
‘none’
- wave, throughput
Wavelength set in user unit and associated unitless emissivity.
- Type:
array_like
- class pysynphot.spectrum.CompositeSpectralElement(component1, component2)
Bases:
SpectralElement
Class to handle composite spectrum involving bandpasses.
- Parameters:
component1 (
SpectralElement
) – Input bandpass.component2 (
SpectralElement
) – Input bandpass.
- component1, component2
Same as inputs.
- isAnalytic
Flag to indicate whether this is an analytic spectrum. This is only
True
if both inputs are analytic.- Type:
- warnings
To store warnings, which are inherited from both input sources. If inputs have the same warning keyword, the one from
component2
is used.- Type:
- primary_area
Area of the telescope. This is inherited from either of the inputs, if available (not
None
). If inputs have different values, an exception is raised.- Type:
number or
None
- binset
This is reserved to be used by
ObsModeBandpass
.- Type:
- waveunits
User unit inherited from inputs, where both inputs are required to have the same unit or an exception will be raised.
- Type:
- wave, throughput
Wavelength set in user unit and associated unitless throughput.
- Type:
array_like
- Raises:
NotImplementedError – Inputs have different wavelength units.
TypeError – Both input spectra must be bandpasses.
pysynphot.exceptions.IncompatibleSources – Input spectra have different telescope areas defined.
- GetWaveSet()
Obtain the wavelength set for the composite spectrum. This is done by using
MergeWaveSets()
to form a union of wavelength sets from its components.- Returns:
waveset – Composite wavelength set.
- Return type:
array_like
- complist()
Return a list of all components and sub-components.
- pysynphot.spectrum.MergeWaveSets(waveset1, waveset2)
Return the union of the two wavelength sets.
The union is computed using
numpy.union1d
, unless one or both of them isNone
.The merged result may sometimes contain numbers which are nearly equal but differ at levels as small as 1E-14. Having values this close together can cause problems due to effectively duplicate wavelength values. Therefore, wavelength values having differences smaller than or equal to
pysynphot.spectrum.MERGETHRESH
(defaults to 1E-12) are considered as the same.
- pysynphot.spectrum.trimSpectrum(sp, minw, maxw)
Create a new spectrum with trimmed upper and lower ranges.
- Parameters:
sp (
SourceSpectrum
) – Spectrum to trim.minw (number) – Lower and upper limits (inclusive) for the wavelength set in the trimmed spectrum.
maxw (number) – Lower and upper limits (inclusive) for the wavelength set in the trimmed spectrum.
- Returns:
result – Trimmed spectrum.
- Return type:
pysynphot.spparser
This module implements the Language Parser. It uses pysynphot.spark.
pysynphot.tables
This module handles graph and component tables. They are discussed in detail in Appendix C: TMG, TMC, and TMT Files.
- class pysynphot.tables.GraphTable(GFile=None)
Class to handle a graph table.
- Parameters:
GFile (str) – Graph table filename.
- keywords
Values from
KEYWORD
column in EXT 1, converted to lowercase.- Type:
array_like
- innodes, outnodes, compnames, thcompnames
Values from
INNODE
,OUTNODE
,COMPNAME
, andTHCOMPNAME
columns in EXT 1.- Type:
array_like
- primary_area
Value from
PRIMAREA
in EXT 0 header, if exists.- Type:
number
- Raises:
TypeError – No filename given.
- GetComponentsFromGT(modes, innode)
Obtain components from graph table for the given observation mode keywords and starting node.
Note
This prints extra information to screen if
pysynphot.tables.DEBUG
is set toTrue
.- Parameters:
- Returns:
components, thcomponents – List of optical and thermal component names.
- Return type:
- Raises:
KeyError – No matches found for one of the keywords.
ValueError – Incomplete observation mode or unused keyword(s) detected.
- class pysynphot.tables.CompTable(CFile=None)
Class to handle a component table.
- Parameters:
CFile (str) – Component table filename.
- name
Same as input
CFile
.
- compnames, filenames
Values from
COMPNAME
andFILENAME
columns in EXT 1.- Type:
array_like
- Raises:
TypeError – No filename given.
pysynphot.units
This module handles wavelength and flux units. Constants used for unit conversion are also defined here.
- pysynphot.units.Units(uname)
Generate a unit object.
- Parameters:
uname (str) – Wavelength or flux unit name.
- Returns:
unit – Unit object.
None
means unitless.- Return type:
- Raises:
ValueError – Unknown unit name.
- class pysynphot.units.BaseUnit(uname)
Base class for all units.
- Parameters:
uname (str) – Unit name.
- name
Same as input
uname
.
- Dispatch
This is used by sub-classes for unit conversion by mapping target unit name to the relevant conversion method.
- Type:
dict or
None
- Convert(wave, flux, target_units)
Perform unit conversion.
Note
This is only applicable to some of the available flux conversions. All other sub-classes must re-implement this method.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
target_units (str) – Unit to convert to.
- Returns:
result – Converted values.
- Return type:
number or array_like
- Raises:
TypeError – Conversion to given unit is not allowed.
- class pysynphot.units.WaveUnits
Bases:
BaseUnit
Base unit for wavelength.
Since Angstrom is the internal unit used by pysynphot calculations, a wavelength unit can always convert to/from Angstrom. Conversion between two arbitrary units uses Angstrom as the intermediate unit.
- Convert(wave, target_units)
Perform unit conversion.
- class pysynphot.units.Angstrom
Bases:
WaveUnits
Class to handle Angstrom unit.
Since Angstrom is the internal wavelength unit, it can convert to all the other supported units.
- ToAngstrom(wave)
Convert to Angstrom.
Since there is no real conversion necessary, this returns a copy of input (if array) or just the input (if scalar). An input array is copied to avoid modifying the input in subsequent pysynphot processing.
- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToCm(wave)
Convert to cm.
\[\mathrm{cm} = 10^{-8} \; \AA\]- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToHz(wave)
Convert to Hz.
\[\mathrm{Hz} = \frac{c}{\AA}\]where \(c\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToInverseMicron(wave)
Convert to inverse micron.
\[\mu \mathrm{m}^{-1} = 10^{4} \; \AA^{-1}\]- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToMeter(wave)
Convert to meter.
\[\mathrm{m} = 10^{-10} \; \AA\]- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToMicron(wave)
Convert to micron.
\[\mu \mathrm{m} = 10^{-4} \; \AA\]- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToMm(wave)
Convert to mm.
\[\mathrm{mm} = 10^{-7} \; \AA\]- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToNm(wave)
Convert to nm.
\[\mathrm{nm} = 10^{-1} \; \AA\]- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.InverseMicron
Bases:
WaveUnits
Class to handle inverse micron unit.
- ToAngstrom(wave)
Convert to Angstrom.
\[\AA = \frac{10^{4}}{\mu \mathrm{m}^{-1}}\]- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units._MetricWavelength
Bases:
WaveUnits
Class to handle meter unit and its prefixes.
Note
Angstrom
is not a sub-class of this because as an internal unit, it requires special handling.- ToAngstrom(wave)
Convert to Angstrom.
Conversion is simply the input values multiplied by a factor specific to its sub-class.
- Parameters:
wave (number or array_like) – Wavelength values to be used for conversion.
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.Nm
Bases:
_MetricWavelength
Class to handle nm unit.
- class pysynphot.units.Micron
Bases:
_MetricWavelength
Class to handle micron unit.
- class pysynphot.units.Mm
Bases:
_MetricWavelength
Class to handle mm unit.
- class pysynphot.units.Cm
Bases:
_MetricWavelength
Class to handle cm unit.
- class pysynphot.units.Meter
Bases:
_MetricWavelength
Class to handle meter unit.
- class pysynphot.units.FluxUnits
Bases:
BaseUnit
Base unit for flux.
Since
photlam
is the internal unit used by pysynphot calculations, a flux unit can always convert to/fromphotlam
. Conversion between two arbitrary units usesphotlam
as the intermediate unit.Note
To support source spectrum renormalization without introducing circular import, all supported flux units must have their
StdSpectrum
attributes defined separately inDefineStdSpectraForUnits()
.- nativewave
Native wavelength unit associated with the flux unit. This is for informational purpose only.
- Type:
- Convert(wave, flux, target_units, area=None)
Perform unit conversion.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
target_units (str) – Unit to convert to.
area (number or
None
) – Telescope area, if applicable. This is only needed for conversions involving counts orobmag
.
- Returns:
result – Converted values.
- Return type:
number or array_like
- Raises:
TypeError – Conversion to given unit is not allowed.
- class pysynphot.units.Photlam
Bases:
FluxUnits
Class to handle
photlam
unit.Since
photlam
is the internal flux unit, it can convert to all the other supported units.\[\mathrm{photlam} = \mathrm{photon} \; \mathrm{s}^{-1} \; \mathrm{cm}^{-2} \; \AA^{-1}\]- ToABMag(wave, flux, **kwargs)
Convert to
abmag
.\[\mathrm{AB}_{\nu} = -2.5 \; \log(h \lambda \; \mathrm{photlam}) - 48.6\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToCounts(wave, flux, area=None)
Convert to counts.
\[\mathrm{counts} = \delta \lambda \; \times \; \mathrm{area} \; \times \; \mathrm{photlam}\]where \(\delta \lambda\) represent bin widths derived from
calculate_bin_edges()
andcalculate_bin_widths()
, using the input wavelength values as bin centers.- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
area (number or
None
) – Telescope collecting area. If not given, default value from Reference Data is used.
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToFlam(wave, flux, **kwargs)
Convert to
flam
.\[\mathrm{flam} = \frac{hc}{\lambda} \; \mathrm{photlam}\]where \(h\) and \(c\) are as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToFnu(wave, flux, **kwargs)
Convert to
fnu
.\[\mathrm{fnu} = h \lambda \; \mathrm{photlam}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToJy(wave, flux, **kwargs)
Convert to Jy.
\[\mathrm{Jy} = 10^{23} h \lambda \; \mathrm{photlam}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToOBMag(wave, flux, area=None)
Convert to
obmag
.\[\mathrm{obmag} = -2.5 \; \log(\delta \lambda \; \times \; \mathrm{area} \; \times \; \mathrm{photlam})\]where \(\delta \lambda\) represent bin widths derived from
calculate_bin_edges()
andcalculate_bin_widths()
, using the input wavelength values as bin centers.- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
area (number or
None
) – Telescope collecting area. If not given, default value from Reference Data is used.
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.Since there is no real conversion necessary, this returns a copy of input flux (if array) or just the input (if scalar). An input array is copied to avoid modifying the input in subsequent pysynphot processing.
- Parameters:
wave (number or array_like) – Wavelength (not used) and flux values to be used for conversion.
flux (number or array_like) – Wavelength (not used) and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToPhotnu(wave, flux, **kwargs)
Convert to
photnu
.\[\mathrm{photnu} = \frac{\lambda^{2}}{c} \; \mathrm{photlam}\]where \(c\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToSTMag(wave, flux, **kwargs)
Convert to
stmag
.\[\mathrm{ST}_{\lambda} = -2.5 \; \log(\frac{hc}{\lambda} \; \mathrm{photlam}) - 21.1\]where \(h\) and \(c\) are as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- ToVegaMag(wave, flux, **kwargs)
Convert to
vegamag
.\[\mathrm{vegamag} = -2.5 \; \log(\frac{\mathrm{photlam}}{f_{\mathrm{Vega}}})\]where \(f_{\mathrm{Vega}}\) is the flux of Vega resampled at given wavelength values and converted to
photlam
.- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- TomJy(wave, flux, **kwargs)
Convert to mJy.
\[\mathrm{mJy} = 10^{26} h \lambda \; \mathrm{photlam}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- TomuJy(wave, flux, **kwargs)
Convert to \(\mu \mathrm{Jy}\).
\[\mu \mathrm{Jy} = 10^{29} h \lambda \; \mathrm{photlam}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- TonJy(wave, flux, **kwargs)
Convert to nJy.
\[\mathrm{nJy} = 10^{32} h \lambda \; \mathrm{photlam}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.Flam
Bases:
FluxUnits
Class to handle
flam
unit.\[\mathrm{flam} = \mathrm{erg} \; \mathrm{s}^{-1} \; \mathrm{cm}^{-2} \; \AA^{-1}\]- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[\mathrm{photlam} = \frac{\lambda}{hc} \; \mathrm{flam}\]where \(h\) and \(c\) are as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.Photnu
Bases:
FluxUnits
Class to handle
photnu
unit.\[\mathrm{photnu} = \mathrm{photon} \; \mathrm{s}^{-1} \; \mathrm{cm}^{-2} \; \mathrm{Hz}^{-1}\]- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[\mathrm{photlam} = \frac{c}{\lambda^{2}} \; \mathrm{photnu}\]where \(c\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.Fnu
Bases:
FluxUnits
Class to handle
fnu
unit.\[\mathrm{fnu} = \mathrm{erg} \; \mathrm{s}^{-1} \mathrm{cm}^{-2} \mathrm{Hz}^{-1}\]- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[\mathrm{photlam} = \frac{1}{h \lambda} \; \mathrm{fnu}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.Jy
Bases:
FluxUnits
Class to handle Jy unit.
\[\mathrm{Jy} = 10^{-23} \; \mathrm{erg} \; \mathrm{s}^{-1} \; \mathrm{cm}^{-2} \mathrm{Hz}^{-1}\]- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[\mathrm{photlam} = \frac{10^{-23}}{h \lambda} \; \mathrm{Jy}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.mJy
Bases:
FluxUnits
Class to handle mJy unit.
\[\mathrm{mJy} = 10^{-26} \; \mathrm{erg} \; \mathrm{s}^{-1} \; \mathrm{cm}^{-2} \mathrm{Hz}^{-1}\]- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[\mathrm{photlam} = \frac{10^{-26}}{h \lambda} \; \mathrm{mJy}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.muJy
Bases:
FluxUnits
Class to handle \(\mu \mathrm{Jy}\) unit.
\[\mu \mathrm{Jy} = 10^{-29} \; \mathrm{erg} \; \mathrm{s}^{-1} \; \mathrm{cm}^{-2} \mathrm{Hz}^{-1}\]- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[\mathrm{photlam} = \frac{10^{-29}}{h \lambda} \; \mu \mathrm{Jy}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.nJy
Bases:
FluxUnits
Class to handle nJy unit.
\[\mathrm{nJy} = 10^{-32} \; \mathrm{erg} \; \mathrm{s}^{-1} \; \mathrm{cm}^{-2} \mathrm{Hz}^{-1}\]- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[\mathrm{photlam} = \frac{10^{-32}}{h \lambda} \; \mathrm{nJy}\]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.Counts
Bases:
FluxUnits
Class to handle counts unit. See Counts and Magnitudes for more details.
- ToPhotlam(wave, flux, area=None)
Convert to
photlam
.\[\mathrm{photlam} = \frac{\mathrm{counts}}{\delta \lambda \; \times \; \mathrm{area}}\]where \(\delta \lambda\) represent bin widths derived from
calculate_bin_edges()
andcalculate_bin_widths()
, using the input wavelength values as bin centers.- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
area (number or
None
) – Telescope collecting area. If not given, default value from Reference Data is used.
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.LogFluxUnits
Bases:
FluxUnits
Base class for magnitude, which often requires special handling.
- name, zeropoint
To be set by sub-class of a specific magnitude unit.
- Type:
- nativewave
Native wavelength unit associated with the flux unit. This is for informational purpose only.
- Type:
- class pysynphot.units.ABMag
Bases:
LogFluxUnits
Class to handle
abmag
unit. See Counts and Magnitudes for more details.- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[ \begin{align}\begin{aligned}m = -0.4 \; (\mathrm{AB}_{\nu} + 48.6)\\\mathrm{photlam} = \frac{10^{m}}{h \lambda}\end{aligned}\end{align} \]where \(h\) is as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.STMag
Bases:
LogFluxUnits
Class to handle
stmag
unit. See Counts and Magnitudes for more details.- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[ \begin{align}\begin{aligned}m = -0.4 \; (\mathrm{ST}_{\lambda} + 21.1)\\\mathrm{photlam} = \frac{10^{m} \lambda}{hc}\end{aligned}\end{align} \]where \(h\) and \(c\) are as defined in Constants.
- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.OBMag
Bases:
LogFluxUnits
Class to handle
obmag
unit. See Counts and Magnitudes for more details.- ToPhotlam(wave, flux, area=None)
Convert to
photlam
.\[\mathrm{photlam} = \frac{10^{-0.4 \; \mathrm{obmag}}}{\delta \lambda \; \times \; \mathrm{area}}\]where \(\delta \lambda\) represent bin widths derived from
calculate_bin_edges()
andcalculate_bin_widths()
, using the input wavelength values as bin centers.- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
area (number or
None
) – Telescope collecting area. If not given, default value from Reference Data is used.
- Returns:
result – Converted values.
- Return type:
number or array_like
- class pysynphot.units.VegaMag
Bases:
LogFluxUnits
Class to handle
vegamag
unit. See Counts and Magnitudes for more details.- linunit, zeropoint
This is not used.
- Type:
- ToPhotlam(wave, flux, **kwargs)
Convert to
photlam
.\[\mathrm{photlam} = 10^{-0.4 \; \mathrm{vegamag}} \; f_{\mathrm{Vega}}\]where \(f_{\mathrm{Vega}}\) is the flux of Vega resampled at given wavelength values and converted to
photlam
.- Parameters:
wave (number or array_like) – Wavelength and flux values to be used for conversion.
flux (number or array_like) – Wavelength and flux values to be used for conversion.
kwargs (dict) – Extra keywords (not used).
- Returns:
result – Converted values.
- Return type:
number or array_like
pysynphot.wavetable
This module handles selection of the appropriate wavelength table for a given observation mode. This is the same selection as used by ETC.
Its wavecat_file
(see below) can contain a mix of the following:
Name of the ASCII file containing the wavelength values. The filename can contain IRAF-style path shortcut. The file must only contain one column, with a single wavelength value in Angstrom in each row. The file can also contain comment lines that begin with “#”, which will be skipped.
A string of comma-separated coefficients that describe how to construct the wavelength table in Angstrom, in the format of
(c0,c1[,c2[,c3]])
, wherec2
andc3
are optional. They are used for the following computation. Basically, the wavelength table runs fromc0
toc1
, with constant \(\delta \lambda\) ifc2
is undefined, or ifc2
is given but notc3
, or variable \(\delta \lambda\) if bothc2
andc3
are given:\[ \begin{align}\begin{aligned}\begin{split}c_{2} = \left \{ \begin{array}{ll} c_{2} & : \mathrm{if given} \\ (c_{1} - c_{0})/1999 & : \mathrm{else, where 1999 was taken from IRAF STSDAS SYNPHOT} \end{array} \right.\end{split}\\\begin{split}c_{3} = \left \{ \begin{array}{ll} c_{3} & : \mathrm{if given} \\ c_{2} & : \mathrm{else} \end{array} \right.\end{split}\\n = \mathrm{int}(\frac{2 \; (c_{1} - c_{0})}{c_{3} + c_{2}}) + 1\\a = \frac{0.25 \; (c_{3}^{2} - c_{2}^{2})}{c_{1} - c_{0}}\\\lambda_{i=0,n-1} = (a i + c_{2}) i + c_{0}\end{aligned}\end{align} \]
Example contents of a wavecat_file
:
# OBSMODE FILENAME_OR_COEFFS
stis,e140h,c1598 (1497.0,1699.0,0.0066,0.0075)
stis,g230l (1568.0,3184.0,1.6)
stis,prism synphot$data/prism.dat
stis,prism,c1200 synphot$data/prism.dat
Global Variables
pysynphot.wavetable.wavecat_file
- This is the same aspysynphot.locations.wavecat
. It is the data file used in this module.pysynphot.wavetable.wavetable
- This is aWavetable
object created usingpysynphot.wavetable.wavecat_file
.
- class pysynphot.wavetable.Wavetable(fname)
Class to handle wavelength table.
__getitem__()
is used to look up the wavelength table. It raisesKeyError
orValueError
if look-up fails. The look-up result is resolved into actual wavelength values bybandWave()
.- Parameters:
fname (str) – Data file.
- file
Same as input
fname
.
- lookup
Look-up table using
obsmode
string as key. This is used by default for direct match.- Type:
- setlookup
Same as
lookup
but theobsmode
string is converted into a frozen set consisting of its components. This is used for partial look-up if there is no direct match.- Type:
- Raises:
ValueError – Failed to parse input file.
Examples
>>> wavetab = S.wavetable.Wavetable(S.wavetable.wavecat_file) >>> wavetab['stis,g230l'] '(1568.0,3184.0,1.6)' >>> wavetab['stis,prism'] 'synphot$data/prism.dat'