Sed¶
- class rubin_sim.phot_utils.Sed(wavelen=None, flambda=None, fnu=None, badval=nan, name=None)¶
Bases:
object
“Hold and use spectral energy distributions (SEDs)
Methods Summary
add_ccm_dust
(a_x, b_x[, a_v, ebv, r_v, ...])Add dust model extinction to the SED, modifying flambda and fnu.
add_dust
(a_x, b_x[, a_v, ebv, r_v, wavelen, ...])Add dust model extinction to the SED, modifying flambda and fnu.
calc_adu
(bandpass, phot_params[, wavelen, fnu])Calculate the number of adu from camera, using sb and fnu.
calc_ergs
(bandpass)Integrate the SED over a bandpass directly.
calc_flux
(bandpass[, wavelen, fnu])Integrate the specific flux density of the object over the normalized response curve of a bandpass, giving a flux in Janskys (10^-23 ergs/s/cm^2/Hz) through the normalized response curve, as detailed in Section 4.1 of the LSST design document LSE-180 and Section 2.6 of the LSST Science Book (http://ww.lsst.org/scientists/scibook).
calc_flux_norm
(magmatch, bandpass[, ...])Calculate the fluxNorm (SED normalization value for a given mag) for a sed.
calc_mag
(bandpass[, wavelen, fnu])Calculate the AB magnitude of an object using the normalized system response (phi from Section 4.1 of the LSST design document LSE-180).
Reset all data in sed to None.
flambda_tofnu
([wavelen, flambda])Convert flambda into fnu.
flux_from_mag
(mag)Convert a magnitude back into a flux (implies knowledge of the zeropoint, which is stored in this class)
fnu_toflambda
([wavelen, fnu])Convert fnu into flambda.
Return copy of wavelen/flambda.
Return copy of wavelen/fnu, without altering self.
mag_from_flux
(flux)Convert a flux into a magnitude (implies knowledge of the zeropoint, which is stored in this class)
many_flux_calc
(phiarray, wavelen_step[, ...])Calculate fluxes of a single sed for which fnu has been evaluated in a set of bandpasses for which phiarray has been set up to have the same wavelength grid as the SED in units of ergs/cm^2/sec.
many_mag_calc
(phiarray, wavelen_step[, ...])Calculate many magnitudes for many bandpasses using a single sed.
multiply_flux_norm
(flux_norm[, wavelen, fnu])Multiply wavelen/fnu (or self.wavelen/fnu) by fluxnorm.
multiply_sed
(other_sed[, wavelen_step])Multiply two SEDs together - flambda * flambda - and return a new sed object.
read_sed_flambda
(filename[, name, cache_sed])Read a file containing [lambda Flambda] (lambda in nm) (Flambda erg/cm^2/s/nm).
read_sed_fnu
(filename[, name])Read a file containing [lambda Fnu] (lambda in nm) (Fnu in Jansky).
redshift_sed
(redshift[, dimming, wavelen, ...])Redshift an SED, optionally adding cosmological dimming.
renormalize_sed
([wavelen, flambda, fnu, ...])Renormalize sed in flambda to have normflux=normvalue @ lambdanorm or averaged over gap.
resample_sed
([wavelen, flux, wavelen_match, ...])Resample flux onto grid defined by min/max/step OR another wavelength array.
set_flat_sed
([wavelen_min, wavelen_max, ...])Populate the wavelength/flambda/fnu fields in sed according to a flat fnu source.
set_sed
(wavelen[, flambda, fnu, name])Populate wavelen/flambda fields in sed by giving lambda/flambda or lambda/fnu array.
setup_cc_mab
([wavelen])Calculate a(x) and b(x) for CCM dust model.
setup_ccm_ab
([wavelen])Calculate a(x) and b(x) for CCM dust model.
setup_o_donnell_ab
([wavelen])Calculate a(x) and b(x) for O'Donnell dust model.
setup_phi_array
(bandpasslist)Sets up a 2-d numpy phi array from bandpasslist suitable for input to Sed's many_mag_calc.
synchronize_sed
([wavelen_min, wavelen_max, ...])Set all wavelen/flambda/fnu values, potentially on min/max/step grid.
write_sed
(filename[, print_header, ...])Write SED (wavelen, flambda, optional fnu) out to file.
Methods Documentation
- add_ccm_dust(a_x, b_x, a_v=None, ebv=None, r_v=3.1, wavelen=None, flambda=None)¶
Add dust model extinction to the SED, modifying flambda and fnu.
Get a_x and b_x either from setupCCMab or setupODonnell_ab
Specify any two of A_V, E(B-V) or R_V (=3.1 default).
- add_dust(a_x, b_x, a_v=None, ebv=None, r_v=3.1, wavelen=None, flambda=None)¶
Add dust model extinction to the SED, modifying flambda and fnu.
Get a_x and b_x either from setupCCMab or setupODonnell_ab
Specify any two of A_V, E(B-V) or R_V (=3.1 default).
- calc_adu(bandpass, phot_params, wavelen=None, fnu=None)¶
Calculate the number of adu from camera, using sb and fnu.
Given wavelen/fnu arrays or use self. Self or passed wavelen/fnu arrays will be unchanged. Calculating the AB mag requires the wavelen/fnu pair to be on the same grid as bandpass; (temporary values of these are used).
- Parameters:
- bandpass
rubin_sim.phot_utils.Bandpass
- phot_params
rubin_sim.phot_utils.PhotometricParameters
- wavelen
np.ndarray
, optional wavelength grid in nm
- fnu
np.ndarray
, optional flux in Janskys
- If wavelen and fnu are not specified, this will just use self.wavelen and
- self.fnu
- bandpass
- calc_ergs(bandpass)¶
Integrate the SED over a bandpass directly. If self.flambda is in ergs/s/cm^2/nm and bandpass.sb is the unitless probability that a photon of a given wavelength will pass through the system, this method will return the ergs/s/cm^2 of the source observed through that bandpass (i.e. it will return the integral
int self.flambda(lambda) * bandpass.sb(lambda) * dlambda
This is to be contrasted with self.calc_flux(), which returns the integral of the source’s specific flux density over the normalized response function of bandpass, giving a flux in Janskys (10^-23 erg/cm^2/s/Hz), which should be though of as a weighted average of the specific flux density of the source over the normalized response function, as detailed in Section 4.1 of the LSST design document LSE-180.
- Parameters:
- bandpass is an instantiation of the Bandpass class
- Returns:
- The flux of the current SED through the bandpass in ergs/s/cm^2
- calc_flux(bandpass, wavelen=None, fnu=None)¶
Integrate the specific flux density of the object over the normalized response curve of a bandpass, giving a flux in Janskys (10^-23 ergs/s/cm^2/Hz) through the normalized response curve, as detailed in Section 4.1 of the LSST design document LSE-180 and Section 2.6 of the LSST Science Book (http://ww.lsst.org/scientists/scibook). This flux in Janskys (which is usually though of as a unit of specific flux density), should be considered a weighted average of the specific flux density over the normalized response curve of the bandpass. Because we are using the normalized response curve (phi in LSE-180), this quantity will depend only on the shape of the response curve, not its absolute normalization.
Note: the way that the normalized response curve has been defined (see equation 5 of LSE-180) is appropriate for photon-counting detectors, not calorimeters.
Passed wavelen/fnu arrays will be unchanged, but if uses self will check if fnu is set.
Calculating the AB mag requires the wavelen/fnu pair to be on the same grid as bandpass; (temporary values of these are used).
- calc_flux_norm(magmatch, bandpass, wavelen=None, fnu=None)¶
Calculate the fluxNorm (SED normalization value for a given mag) for a sed.
Equivalent to adjusting a particular f_nu to Jansky’s appropriate for the desired mag. Can pass wavelen/fnu or apply to self.
- calc_mag(bandpass, wavelen=None, fnu=None)¶
Calculate the AB magnitude of an object using the normalized system response (phi from Section 4.1 of the LSST design document LSE-180).
Can pass wavelen/fnu arrays or use self. Self or passed wavelen/fnu arrays will be unchanged. Calculating the AB mag requires the wavelen/fnu pair to be on the same grid as bandpass; (but only temporary values of these are used).
- clear_sed()¶
Reset all data in sed to None.
- flambda_tofnu(wavelen=None, flambda=None)¶
Convert flambda into fnu.
This routine assumes that flambda is in ergs/cm^s/s/nm and produces fnu in Jansky. Can act on self or user can provide wavelen/flambda and get back wavelen/fnu.
- flux_from_mag(mag)¶
Convert a magnitude back into a flux (implies knowledge of the zeropoint, which is stored in this class)
- fnu_toflambda(wavelen=None, fnu=None)¶
Convert fnu into flambda.
Assumes fnu in units of Jansky and flambda in ergs/cm^s/s/nm. Can act on self or user can give wavelen/fnu and get wavelen/flambda returned.
- get_sed_flambda()¶
Return copy of wavelen/flambda.
- get_sed_fnu()¶
Return copy of wavelen/fnu, without altering self.
- mag_from_flux(flux)¶
Convert a flux into a magnitude (implies knowledge of the zeropoint, which is stored in this class)
- many_flux_calc(phiarray, wavelen_step, observed_bandpass_ind=None)¶
Calculate fluxes of a single sed for which fnu has been evaluated in a set of bandpasses for which phiarray has been set up to have the same wavelength grid as the SED in units of ergs/cm^2/sec. It is assumed that
self.fnu
is set before calling this method, and that phiArray has the same wavelength grid as the Sed.- Parameters:
- phiarray: `np.ndarray`, mandatory
phiarray corresponding to the list of bandpasses in which the band fluxes need to be calculated, in the same wavelength grid as the SED
- wavelen_step: `float`, mandatory
the uniform grid size of the SED
- observed_bandpass_ind: list of integers, optional, defaults to None
list of indices of phiarray corresponding to observed bandpasses, if None, the original phiarray is returned
- Returns:
np.ndarray
with size equal to number of bandpass filters band flux- values in units of ergs/cm^2/sec
- grid as the Sed and that
sed.fnu
has been calculated for the sed, - perhaps using
sed.flambda_tofnu()
. This requires calling sed.setupPhiArray()
first. These assumptions are to avoid error- checking within this function (for speed), but could lead to errors if
- method is used incorrectly.
- Note on units: Fluxes calculated this way will be the flux density integrated over the
- weighted response curve of the bandpass. See equaiton 2.1 of the LSST Science Book
- http://www.lsst.org/scientists/scibook
- many_mag_calc(phiarray, wavelen_step, observed_bandpass_ind=None)¶
Calculate many magnitudes for many bandpasses using a single sed.
This method assumes that there will be flux within a particular bandpass (could return ‘-Inf’ for a magnitude if there is none). Use setupPhiArray first, and note that Sed.many_mag_calc assumes phiArray has the same wavelength grid as the Sed, and that fnu has already been calculated for Sed. These assumptions are to avoid error checking within this function (for speed), but could lead to errors if method is used incorrectly.
- Parameters:
- phiarray: `np.ndarray`, mandatory
phiarray corresponding to the list of bandpasses in which the band fluxes need to be calculated, in the same wavelength grid as the SED
- wavelen_step: `float`, mandatory
the uniform grid size of the SED
- observed_bandpass_ind: list of integers, optional, defaults to None
list of indices of phiarray corresponding to observed bandpasses, if None, the original phiarray is returned
- multiply_flux_norm(flux_norm, wavelen=None, fnu=None)¶
Multiply wavelen/fnu (or self.wavelen/fnu) by fluxnorm.
Returns wavelen/fnu arrays (or updates self). Note that multiply_flux_norm does not regrid self.wavelen/flambda/fnu at all.
- multiply_sed(other_sed, wavelen_step=None)¶
Multiply two SEDs together - flambda * flambda - and return a new sed object.
Unless the two wavelength arrays are equal, returns a SED gridded with stepsize wavelen_step over intersecting wavelength region. Does not alter self or other_sed.
- read_sed_flambda(filename, name=None, cache_sed=True)¶
Read a file containing [lambda Flambda] (lambda in nm) (Flambda erg/cm^2/s/nm).
Does not resample wavelen/flambda onto grid; leave fnu=None.
- read_sed_fnu(filename, name=None)¶
Read a file containing [lambda Fnu] (lambda in nm) (Fnu in Jansky).
Does not resample wavelen/fnu/flambda onto a grid; leaves fnu set.
- redshift_sed(redshift, dimming=False, wavelen=None, flambda=None)¶
Redshift an SED, optionally adding cosmological dimming.
Pass wavelen/flambda or redshift/update self.wavelen/flambda (unsets fnu).
- renormalize_sed(wavelen=None, flambda=None, fnu=None, lambdanorm=500, normvalue=1, gap=0, normflux='flambda', wavelen_step=None)¶
Renormalize sed in flambda to have normflux=normvalue @ lambdanorm or averaged over gap.
Can normalized in flambda or fnu values. wavelen_step specifies the wavelength spacing when using ‘gap’.
Either returns wavelen/flambda values or updates self.
- resample_sed(wavelen=None, flux=None, wavelen_match=None, wavelen_min=None, wavelen_max=None, wavelen_step=None, force=False)¶
Resample flux onto grid defined by min/max/step OR another wavelength array.
Give method wavelen/flux OR default to self.wavelen/self.flambda. Method either returns wavelen/flambda (if given those arrays) or updates wavelen/flambda in self. If updating self, resets fnu to None. Method will first check if resampling needs to be done or not, unless ‘force’ is True.
- set_flat_sed(wavelen_min=None, wavelen_max=None, wavelen_step=None, name='Flat')¶
Populate the wavelength/flambda/fnu fields in sed according to a flat fnu source.
- set_sed(wavelen, flambda=None, fnu=None, name='FromArray')¶
Populate wavelen/flambda fields in sed by giving lambda/flambda or lambda/fnu array.
If flambda present, this overrides fnu. Method sets fnu=None unless only fnu is given. Sets wavelen/flambda or wavelen/flambda/fnu over wavelength array given.
- setup_cc_mab(wavelen=None)¶
Calculate a(x) and b(x) for CCM dust model. (x=1/wavelen).
If wavelen not specified, calculates a and b on the own object’s wavelength grid. Returns a(x) and b(x) can be common to many seds, wavelen is the same.
This method sets up extinction due to the model of Cardelli, Clayton and Mathis 1989 (ApJ 345, 245)
- setup_ccm_ab(wavelen=None)¶
Calculate a(x) and b(x) for CCM dust model. (x=1/wavelen).
If wavelen not specified, calculates a and b on the own object’s wavelength grid. Returns a(x) and b(x) can be common to many seds, wavelen is the same.
This method sets up extinction due to the model of Cardelli, Clayton and Mathis 1989 (ApJ 345, 245)
- setup_o_donnell_ab(wavelen=None)¶
Calculate a(x) and b(x) for O’Donnell dust model. (x=1/wavelen).
If wavelen not specified, calculates a and b on the own object’s wavelength grid. Returns a(x) and b(x) can be common to many seds, wavelen is the same.
This method sets up the extinction parameters from the model of O’Donnel 1994 (ApJ 422, 158)
- setup_phi_array(bandpasslist)¶
Sets up a 2-d numpy phi array from bandpasslist suitable for input to Sed’s many_mag_calc.
This is intended to be used once, most likely before using Sed’s many_mag_calc many times on many SEDs. Returns 2-d phi array and the wavelen_step (dlambda) appropriate for that array.
- synchronize_sed(wavelen_min=None, wavelen_max=None, wavelen_step=None)¶
Set all wavelen/flambda/fnu values, potentially on min/max/step grid.
Uses flambda to recalculate fnu. If wavelen min/max/step are given, resamples wavelength/flambda/fnu onto an even grid with these values.
- write_sed(filename, print_header=None, print_fnu=False, wavelen_min=None, wavelen_max=None, wavelen_step=None)¶
Write SED (wavelen, flambda, optional fnu) out to file.
Option of adding a header line (such as version info) to output file. Does not alter self, regardless of grid or presence/absence of fnu.