pypeit.core.flux_calib module

Module for fluxing routines

pypeit.core.flux_calib.Flam_to_Nlam(wave, zeropoint, zp_min=5.0, zp_max=30.0)[source]

The factor that when multiplied into F_lam converts to N_lam, i.e. 1/S_lam where S_lam equiv F_lam/N_lam

Parameters
  • wave (numpy.ndarray) – Wavelength array, float, shape (nspec,)

  • zeropoint (numpy.ndarray) – zeropoint array, float, shape (nspec,)

Returns

factor – Factor that when multiplied into F_lam converts to N_lam

Return type

numpy.ndarray

pypeit.core.flux_calib.Nlam_to_Flam(wave, zeropoint, zp_min=5.0, zp_max=30.0)[source]

The factor that when multiplied into N_lam converts to F_lam, i.e. S_lam where S_lam equiv F_lam/N_lam

Parameters
  • wave (numpy.ndarray) – Wavelength vector for zeropoint

  • zeropoint (numpy.ndarray) – zeropoint

  • zp_min (float, optional) – Minimum allowed value of the ZP. For smaller values the S_lam factor is set to zero

  • zp_max (float, optional) – Maximum allowed value of the ZP. For larger values the S_lam factor is set to zero

Returns

factor – S_lam factor

Return type

numpy.ndarray

pypeit.core.flux_calib.blackbody_func(a, teff)[source]

Generate a blackbody spectrum based on the normalisation and effective temperature. See Suzuki & Fukugita, 2018, AJ, 156, 219: https://ui.adsabs.harvard.edu/abs/2018AJ….156..219S/abstract

Parameters
  • a (float) – flux normalisation factor (dimensionless)

  • teff (float) – Effective temperature of the blackbody (in units of K)

Returns

  • waves (numpy.ndarray) – wavelengths

  • flam (numpy.ndarray) – flux in units of erg/s/cm^2/A

pypeit.core.flux_calib.compute_zeropoint(wave, N_lam, N_lam_gpm, flam_std_star, tellmodel=None)[source]

Routine to compute the zeropoint and zeropoint_gpm from the N_lam (counts/s/A) of a standard star

Parameters
  • wave (numpy.ndarray) – Wavelength array, float, shape (nspec,)

  • N_lam (numpy.ndarray) – N_lam spectrum of standard star, float, shape (nspec,)

  • N_lam_gpm (numpy.ndarray) – N_lam mask, good pixel mask, boolean, shape (nspec,)

  • flam_std_star (numpy.ndarray) – True standard star spectrum units set of PYPEIT_FLUX_SCALE erg/s/cm^2/sm/Angstrom

  • tellmodel (numpy.ndarray) – Telluric absorption model, optional, shape (nspec,)

Returns

  • zeropoint (numpy.ndarray) – Spectroscopic zeropoint, float, shape (nspec,)

  • zeropoint_gpm (numpy.ndarray) – Zeropoint good pixel mask, bool, shape (nspec,)

pypeit.core.flux_calib.counts2Nlam(wave, counts, counts_ivar, counts_mask, exptime, airmass, longitude, latitude, extinctfilepar)[source]

Convert counts to counts/s/Angstrom Used for flux calibration and to apply extinction correction

Parameters
  • wave (numpy.ndarray) – Wavelength of the star. Shape (nspec,)

  • counts (numpy.ndarray) – Flux (in counts) of the star. Shape (nspec,)

  • counts_ivar (numpy.ndarray) – Inverse variance of the star counts. Shape (nspec,)

  • counts_mask (numpy.ndarray) – Good pixel mask for the counts.

  • exptime (float) – Exposure time in seconds

  • airmass (float) – Airmass

  • longitude (float) – Telescope longitude, used for extinction correction.

  • latitude (float) – Telescope latitude, used for extinction correction

  • extinctfilepar (str) – [sensfunc][UVIS][extinct_file] parameter Used for extinction correction

Returns

Three items:

Return type

tuple

pypeit.core.flux_calib.extinction_correction(wave, airmass, extinct)[source]

Derive extinction correction Based on algorithm in LowRedux (long_extinct)

Parameters
  • wave (numpy.ndarray) – Wavelengths for interpolation. Should be sorted. Assumes angstroms.

  • airmass (float) – Airmass

  • extinct (Table) – Table of extinction values

Returns

flux_corr – Multiplucative flux correction factors at the input wavelengths. i.e. true_flux = correction_factor*observed_flux

Return type

numpy.ndarray

pypeit.core.flux_calib.find_standard(specobj_list)[source]

Routine to identify the standard star given a list of spectra

Take the median boxcar and then take the max flux object (in BOX_COUNTS) as the standard

Parameters

specobj_list (list) – pypeit.specobj.SpecObj list

Returns

mxix – Index of the standard star in the list

Return type

int

pypeit.core.flux_calib.find_standard_file(ra, dec, toler=<Quantity 20. arcmin>, check=False)[source]

Find a match for the input file to one of the archived standard star files (hopefully).

Priority is set by the following search order:

['xshooter', 'calspec', 'esofil', 'noao', 'ing', 'blackbody']
Parameters
  • ra (float) – Object right-ascension in decimal deg

  • dec (float) – Object declination in decimal deg

  • toler (astropy.units.quantity.Quantity, optional) – Tolerance on matching archived standards to input. Expected to be in arcmin.

  • check (bool, optional) – If True, the routine will only check to see if a standard star exists within the input ra, dec, and toler range.

Returns

star_dict – If check is True, return True or False depending on if the object is matched to a library standard star. If check is False and no match is found, return None. Otherwise, return a dictionary with the matching standard star with the following meta data:

  • cal_file: str – Filename table

  • name: str – Star name

  • std_ra: float – RA(J2000)

  • std_dec: float – DEC(J2000)

Return type

dict, bool or None

pypeit.core.flux_calib.fit_zeropoint(wave, Nlam_star, Nlam_ivar_star, gpm_star, std_dict, mask_hydrogen_lines=True, mask_helium_lines=False, polyorder=4, hydrogen_mask_wid=10.0, nresln=20.0, resolution=3000.0, trans_thresh=0.9, polycorrect=True, polyfunc=False, debug=False)[source]

Function to generate the sensitivity function. This function fits a bspline to the 2.5*log10(flux_std/flux_counts). The break points spacing, which determines the scale of variation of the sensitivity function is determined by the nresln parameter.

Parameters
  • wave (numpy.ndarray) – Wavelength of the star. Shape (nspec,)

  • Nlam_star (numpy.ndarray) – counts/second/Angstrom

  • Nlam_ivar_star (numpy.ndarray) – Inverse variance of Nlam_star

  • gpm_star (numpy.ndarray) – Good pixel mask for Nlam_star

  • std_dict (dict) – Dictionary containing information about the standard star returned by flux_calib.get_standard_spectrum

  • mask_hydrogen_lines (bool, optional) – If True, mask stellar hydrogen absorption lines before fitting sensitivity function. Default = True

  • mask_helium_lines (bool, optional) – If True, mask stellar helium absorption lines before fitting sensitivity function. Default = False

  • hydrogen_mask_wid (float, optional) – Parameter describing the width of the mask for or stellar absorption lines (i.e., mask_hydrogen_lines=True) in Angstroms. A region equal to hydrogen_mask_wid on either side of the line center is masked. Default = 10A

  • polycorrect (bool, optional) – Whether you want to interpolate the zeropoint with polynomial in the stellar absortion line regions before fitting with the bspline

  • nresln (float, optional) – Parameter governing the spacing of the bspline breakpoints. default = 20.0

  • resolution (float, optional) – Expected resolution of the standard star spectrum. This should probably be determined from the grating, but is currently hard wired. default=3000.0

  • trans_thresh (float, optional) – Parameter for selecting telluric regions which are masked. Locations below this transmission value are masked. If you have significant telluric absorption you should be using telluric.sensnfunc_telluric. default = 0.9

  • polyfunc (bool, optional) – If True, the zeropoint was a polynomial and not a bspline

Returns

  • zeropoint_data (numpy.ndarray) – Sensitivity function with same shape as wave (nspec,)

  • zeropoint_data_gpm (numpy.ndarray) – Good pixel mask for sensitivity function with same shape as wave (nspec,)

  • zeropoint_fit (numpy.ndarray) – Fitted sensitivity function with same shape as wave (nspec,)

  • zeropoint_fit_gpm (numpy.ndarray) – Good pixel mask for fitted sensitivity function with same shape as wave (nspec,)

Return type

tuple

pypeit.core.flux_calib.get_mask(wave_star, flux_star, ivar_star, mask_star, mask_hydrogen_lines=True, mask_helium_lines=False, mask_telluric=True, hydrogen_mask_wid=10.0, trans_thresh=0.9)[source]

Generate a set of masks from your observed standard spectrum. e.g. Balmer absorption

Parameters
  • wave_star (numpy.ndarray) – wavelength array of your spectrum

  • flux_star (numpy.ndarray) – flux array of your spectrum

  • ivar_star (numpy.ndarray) – ivar array of your spectrum

  • mask_star (bool, optional) – whether you need to mask Hydrogen recombination line region. If False, the returned msk_star are all good.

  • mask_hydrogen_lines (bool, optional) – whether you need to mask hydrogen absorption lines, mask width set by hydrogen_mask_wid

  • mask_helium_lines (bool, optional) – whether you need to mask hydrogen absorption lines, mask width set to \(0.5 \times\) hydrogen_mask_wid

  • mask_telluric (bool, optional) – whether you need to mask telluric region. If False, the returned msk_tell are all good.

  • hydrogen_mask_wid (float, optional) – in units of angstrom Mask parameter for hydrogen recombination absorption lines. A region equal to hydrogen_mask_wid on either side of the line center is masked.

  • trans_thresh (float, optional) – parameter for selecting telluric regions.

Returns

  • msk_bad (bool numpy.ndarray) – mask for bad pixels.

  • msk_star (bool numpy.ndarray) – mask for recombination lines in star spectrum.

  • msk_tell (bool numpy.ndarray) – mask for telluric regions.

pypeit.core.flux_calib.get_sensfunc_factor(wave, wave_zp, zeropoint, exptime, tellmodel=None, extinct_correct=False, airmass=None, longitude=None, latitude=None, extinctfilepar=None, extrap_sens=False)[source]

Get the final sensitivity function factor that will be multiplied into a spectrum in units of counts to flux calibrate it. This code interpolates the sensitivity function and can also multiply in extinction and telluric corrections.

FLAM, FLAM_SIG, and FLAM_IVAR are generated

Parameters
  • wave (float numpy.ndarray) – shape = (nspec,) Senstivity

  • wave_zp (float numpy.ndarray) – Zerooint wavelength vector shape = (nsens,)

  • zeropoint (float numpy.ndarray) – shape = (nsens,) Zeropoint, i.e. sensitivity function

  • exptime (float) –

  • tellmodel (float numpy.ndarray, optional) – shape = (nspec,) Apply telluric correction if it is passed it. Note this is deprecated.

  • extinct_correct (bool, optional) – If True perform an extinction correction. Deafult = False

  • airmass (float, optional) – Airmass used if extinct_correct=True. This is required if extinct_correct=True

  • longitude (float, optional) – longitude in degree for observatory Required for extinction correction

  • latitude – latitude in degree for observatory Required for extinction correction

  • extinctfilepar (str) – [sensfunc][UVIS][extinct_file] parameter Used for extinction correction

  • extrap_sens (bool, optional) – Extrapolate the sensitivity function (instead of crashing out)

Returns

sensfunc_factor – This quantity is defined to be sensfunc_interp/exptime/delta_wave. shape = (nspec,)

Return type

numpy.ndarray

pypeit.core.flux_calib.get_standard_spectrum(star_type=None, star_mag=None, ra=None, dec=None)[source]

Get the standard spetrum using given information of your standard/telluric star.

Parameters
  • star_type (str, optional) – Spectral type of your standard/telluric star

  • star_mag (float, optional) – Apparent magnitude of the telluric star

  • ra (float or str, optional) – Standard right-ascension in decimal degrees (float) -OR- Standard right-ascension in hh:mm:ss string format (e.g.,’05:06:36.6’).

  • dec (float or str, optional) – Standard declination in decimal degrees (float) -OR- Standard declination in dd:mm:ss string format (e.g., 52:52:01.0’)

Returns

std_dict – Dictionary containing the information you provided and the standard/telluric spectrum.

Return type

dict

pypeit.core.flux_calib.load_extinction_data(longitude, latitude, extinctfilepar, toler=<Quantity 5. deg>)[source]

Find the best extinction file to use, based on longitude and latitude. Loads it and returns a Table

Parameters
  • longitude (Geocentric coordinates in degrees (floats).) –

  • latitude (Geocentric coordinates in degrees (floats).) –

  • extinctfilepar ((str)) – The sensfunc[‘extinct_file’] parameter, used to determine which extinction file to load.

  • toler (Angle, optional) – Tolerance for matching detector to site (5 deg)

Returns

ext_file – astropy Table containing the ‘wavelength’, ‘extinct’ data for AM=1.

Return type

astropy.table.Table

pypeit.core.flux_calib.load_filter_file(filter)[source]

Load a system response curve for a given filter. All supported filters can be found at pypeit.data.filters

Parameters

filter (str) – Name of filter

Returns

pypeit.core.flux_calib.mAB_to_cgs(mAB, wvl)[source]

Convert AB magnitudes to flambda cgs unit erg cm^-2 s^-1 A^-1

Parameters
Returns

flux density – f_lambda flux in cgs units

Return type

float or numpy.ndarray

pypeit.core.flux_calib.mask_stellar_helium(wave_star, mask_width=5.0, mask_star=None)[source]

Routine to mask stellar helium recombination lines

Note

This function is pulled out separate from get_mask() because it is used in the telluric module, independent of the remainder of the functionality in get_mask().

Parameters
  • wave_star (numpy.ndarray) – Wavelength of the stellar spectrum shape (nspec,) or (nspec, nimgs)

  • mask_width (float, optional) – width to mask on either side of each line center in Angstroms

  • mask_star (numpy.ndarray, optional) – Incoming star mask to which to add the ionized helium lines (Default: None)

Returns

boolean mask. Same shape as wave_star, True=Good (i.e. does not hit a stellar absorption line).

Return type

numpy.ndarray

pypeit.core.flux_calib.mask_stellar_hydrogen(wave_star, mask_width=10.0, mask_star=None)[source]

Routine to mask stellar hydrogen recombination lines

Note

This function is pulled out separate from get_mask() because it is used in the telluric module, independent of the remainder of the functionality in get_mask().

Parameters
  • wave_star (numpy.ndarray) – Wavelength of the stellar spectrum shape (nspec,) or (nspec, nimgs)

  • mask_width (float, optional) – width to mask on either side of each line center in Angstroms

  • mask_star (numpy.ndarray, optional) – Incoming star mask to which to add the hydrogen lines (Default: None)

Returns

boolean mask. Same shape as wave_star, True=Good (i.e. does not hit a stellar absorption line).

Return type

numpy.ndarray

pypeit.core.flux_calib.scale_in_filter(wave, flux, gpm, scale_dict)[source]

Scale spectra to input magnitude in a given filter

Parameters
Returns

scale – scale value for the flux, i.e. newflux = flux * scale

Return type

float

pypeit.core.flux_calib.sensfunc(wave, counts, counts_ivar, counts_mask, exptime, airmass, std_dict, longitude, latitude, extinctfilepar, ech_orders=None, mask_hydrogen_lines=True, mask_helium_lines=False, polyorder=4, hydrogen_mask_wid=10.0, nresln=20.0, resolution=3000.0, trans_thresh=0.9, polycorrect=True, polyfunc=False, debug=False)[source]

Function to generate the sensitivity function. This function fits a bspline to the 2.5*log10(flux_std/flux_counts). The break points spacing, which determines the scale of variation of the sensitivity function is determined by the nresln parameter. This code can work in different regimes, but NOTE THAT TELLURIC MODE IS DEPRECATED, use telluric.sensfunc_telluric instead

Parameters
  • wave (numpy.ndarray) – Wavelength of the star. Shape (nspec,) or (nspec, norders)

  • counts (ndarray) – Flux (in counts) of the star. Shape (nspec,) or (nspec, norders)

  • counts_ivar (numpy.ndarray) – Inverse variance of the star counts. Shape (nspec,) or (nspec, norders)

  • counts_mask (numpy.ndarray) – Good pixel mask for the counts. Shape (nspec,) or (nspec, norders)

  • exptime (float) – Exposure time in seconds

  • airmass (float) – Airmass

  • std_dict (dict) – Dictionary containing information about the standard star returned by flux_calib.get_standard_spectrum

  • longitude (float) – Telescope longitude, used for extinction correction.

  • latitude (float) – Telescope latitude, used for extinction correction

  • extinctfilepar (str) – [sensfunc][UVIS][extinct_file] parameter Used for extinction correction

  • ech_orders (int numpy.ndarray) – If passed the echelle orders will be added to the meta_table. ech_orders must be a numpy array of integers with the shape (norders,) giving the order numbers

  • mask_hydrogen_lines (bool) – If True, mask stellar hydrogen absorption lines before fitting sensitivity function. Default = True

  • mask_helium_lines (bool) – If True, mask stellar helium absorption lines before fitting sensitivity function. Default = False

  • balm_mask_wid (float) – Parameter describing the width of the mask for or stellar absorption lines (i.e. mask_hydrogen_lines=True). A region equal to balm_mask_wid*resln is masked where resln is the estimate for the spectral resolution in pixels per resolution element.

  • polycorrect (bool) – Whether you want to interpolate the sensfunc with polynomial in the stellar absortion line regions before fitting with the bspline

  • nresln (float) – Parameter governing the spacing of the bspline breakpoints. default = 20.0

  • resolution (float) – Expected resolution of the standard star spectrum. This should probably be determined from the grating, but is currently hard wired. default=3000.0

  • trans_thresh (float) – Parameter for selecting telluric regions which are masked. Locations below this transmission value are masked. If you have significant telluric absorption you should be using telluric.sensnfunc_telluric. default = 0.9

Returns

Returns:

Return type

Tuple

Returns

pypeit.core.flux_calib.standard_zeropoint(wave, Nlam, Nlam_ivar, Nlam_gpm, flam_true, mask_recomb=None, mask_tell=None, maxiter=35, upper=3.0, lower=3.0, func='polynomial', polyorder=5, balm_mask_wid=50.0, nresln=20.0, resolution=2700.0, polycorrect=True, debug=False, polyfunc=False)[source]

Generate a sensitivity function based on observed flux and standard spectrum.

Parameters
  • wave (numpy.ndarray) – wavelength as observed

  • Nlam (numpy.ndarray) – counts/s/Angstrom as observed

  • Nlam_ivar (numpy.ndarray) – inverse variance of counts/s/Angstrom

  • Nlam_gpm (numpy.ndarray) – mask for bad pixels. True is good.

  • flam_true (Quantity array) – standard star true flux (erg/s/cm^2/A)

  • mask_recomb (numpy.ndarray) – mask for hydrogen (and/or helium II) recombination lines. True is good.

  • mask_tell (numpy.ndarray) – mask for telluric regions. True is good.

  • maxiter (integer) – maximum number of iterations for polynomial fit

  • upper (integer) – number of sigma for rejection in polynomial

  • lower (integer) – number of sigma for rejection in polynomial

  • polyorder (integer) – order of polynomial fit

  • balm_mask_wid (float) – Mask parameter for Balmer absorption. A region equal to balm_mask_wid in units of angstrom is masked.

  • nresln (integer/float) – number of resolution elements between breakpoints

  • resolution (integer/float.) – The spectral resolution. This paramters should be removed in the future. The resolution should be estimated from spectra directly.

  • debug (bool) – if True shows some dubugging plots

Returns

  • zeropoint_data (numpy.ndarray) – Sensitivity function with same shape as wave (nspec,)

  • zeropoint_data_gpm (numpy.ndarray) – Good pixel mask for sensitivity function with same shape as wave (nspec,)

  • zeropoint_fit (numpy.ndarray) – Fitted sensitivity function with same shape as wave (nspec,)

  • zeropoint_fit_gpm (numpy.ndarray) – Good pixel mask for fitted sensitivity function with same shape as wave (nspec,)

pypeit.core.flux_calib.stellar_model(V, sptype)[source]

Get the Kurucz stellar model for a given apparent magnitude and spectral type of your standard star. This routine first get the temperature, logg, and bolometric luminosity from the Schmidt-Kaler (1982) table for the given spectral type. It then find the nearest neighbour in the Kurucz stellar atmosphere ATLAS. Finally, the wavelength was converted to Angstrom and the flux density (cgs units) was calculated.

Parameters
  • V (float) – Apparent magnitude of the standard star

  • sptype (str) – Spectral type of the standard star

Returns

pypeit.core.flux_calib.zeropoint_qa_plot(wave, zeropoint_data, zeropoint_data_gpm, zeropoint_fit, zeropoint_fit_gpm, title='Zeropoint QA', axis=None, show=False)[source]

QA plot for zeropoint

Parameters
  • wave (numpy.ndarray) – Wavelength array

  • zeropoint_data (numpy.ndarray) – Zeropoint data array

  • zeropoint_data_gpm (bool numpy.ndarray) – Good pixel mask array for zeropoint_data

  • zeropoint_fit (numpy.ndarray) – Zeropoint fitting array

  • zeropoint_fit_gpm (bool zeropoint_fit) – Good pixel mask array for zeropoint_fit

  • title (str, optional) – Title for the QA plot

  • axis (None or matplotlib.pyplot axis, optional) – axis used for ploting.

  • show (bool, optional) – Whether to show the QA plot

pypeit.core.flux_calib.zeropoint_to_throughput(wave, zeropoint, eff_aperture)[source]

Routine to compute the spectrograph throughput from the zeropoint and effective aperture.

Parameters
  • wave (numpy.ndarray) – Wavelength array shape (nspec,) or (nspec, norders)

  • zeropoint (numpy.ndarray) – Zeropoint array shape (nspec,) or (nspec, norders)

  • eff_aperture (float) – Effective aperture of the telescope in m^2. See spectrograph object

Returns

throughput – Throughput of the spectroscopic setup. Same shape as wave and zeropoint

Return type

numpy.ndarray

pypeit.core.flux_calib.zp_unit_const()[source]

This constant defines the units for the spectroscopic zeropoint. See Flux Calibration.