Skybrightness API

class rubin_sim.skybrightness.Airglow(comp_name='Airglow', sorted_order=['airmass', 'solarFlux'], mags=False)[source]

Bases: BaseSingleInterp

Interpolate the spectra caused by airglow.

class rubin_sim.skybrightness.BaseSingleInterp(comp_name=None, sorted_order=['airmass', 'nightTimes'], mags=False)[source]

Bases: object

Base class for interpolating sky components which only depend on airmass.

Parameters:
  • comp_name (str, optional) – Component name.

  • sorted_order (list [str], optional) – Order of the dimensions in the input .npz files.

  • mags (bool, optional) – Return magnitudes (only) rather than the full spectrum.

__call__(interp_points, filter_names=['u', 'g', 'r', 'i', 'z', 'y'])[source]

At `interp_points (e.g. airmass), return values.

class rubin_sim.skybrightness.LowerAtm(comp_name='LowerAtm', mags=False)[source]

Bases: BaseSingleInterp

Interpolate the spectra caused by the lower atmosphere.

class rubin_sim.skybrightness.MergedSpec(comp_name='MergedSpec', mags=False)[source]

Bases: BaseSingleInterp

Interpolate the combined spectra caused by the sum of the scattered starlight, air glow, upper and lower atmosphere.

class rubin_sim.skybrightness.MoonInterp(comp_name='Moon', sorted_order=['moonSunSep', 'moonAltitude', 'hpid'], mags=False)[source]

Bases: BaseSingleInterp

Read in the saved Lunar spectra and interpolate.

class rubin_sim.skybrightness.ScatteredStar(comp_name='ScatteredStarLight', mags=False)[source]

Bases: BaseSingleInterp

Interpolate the spectra caused by scattered starlight.

class rubin_sim.skybrightness.TwilightInterp(mags=False, dark_sky_mags=None, fit_results=None)[source]

Bases: object

Use the Solar Spectrum to provide an interpolated spectra or magnitudes for the twilight sky.

Parameters:
  • mags (bool) – If True, only return the LSST filter magnitudes, otherwise return the full spectrum

  • dark_sky_mags (dict) – Dict of the zenith dark sky values to be assumed. The twilight fits are done relative to the dark sky level.

  • fit_results (dict) – Dict of twilight parameters based on twilight_func. Keys should be filter names.

interp_mag(interp_points, max_am=3.0, limits=(0.2617993877991494, -0.3490658503988659), filter_names=['u', 'g', 'r', 'i', 'z', 'y'])[source]
Parameters:
  • interp_points (np.ndarray, (N, 3)) – Interpolation points. Should contain sunAlt, airmass and azRelSun.

  • max_am (float, optional) – Maximum airmass to calculate twilight sky to.

  • limits (np.ndarray, (N,), optional) – Sun altitude limits

Returns:

spectra, wavelength

Return type:

np.ndarray, (N, 3), np.ndarray, (M,)

Note

Originally fit the twilight with a cutoff of sun altitude of -11 degrees. I think it can be safely extrapolated farther, but be warned you may be entering a regime where it breaks down.

interp_spec(interp_points, max_am=3.0, limits=(0.2617993877991494, -0.3490658503988659))[source]
Parameters:
  • interp_points (np.ndarray, (N, 3)) – Interpolation points. Should contain sunAlt, airmass and azRelSun.

  • max_am (float, optional) – Maximum airmass to calculate twilight sky to.

  • limits (np.ndarray, (N,), optional) – Sun altitude limits

Returns:

spectra, wavelength

Return type:

np.ndarray, (N, 3), np.ndarray, (M,)

Note

Originally fit the twilight with a cutoff of sun altitude of -11 degrees. I think it can be safely extrapolated farther, but be warned you may be entering a regime where it breaks down.

print_fits_used()[source]

Print out the fit parameters being used

class rubin_sim.skybrightness.UpperAtm(comp_name='UpperAtm', mags=False)[source]

Bases: BaseSingleInterp

Interpolate the spectra caused by the upper atmosphere.

class rubin_sim.skybrightness.ZodiacalInterp(comp_name='Zodiacal', sorted_order=['airmass', 'hpid'], mags=False)[source]

Bases: BaseSingleInterp

Interpolate the zodiacal light based on the airmass and the healpix ID where the healpixels are in ecliptic coordinates, with the sun at ecliptic longitude zero

rubin_sim.skybrightness.all_sky_db(date_id, sql_q=None, dtypes=None, db_address=None, filt='R')[source]

Fetch star and sky magnitudes from a processed all-sky sqlite database.

Parameters:
  • date_id (float) – Date (MJD) to fetch star observation information from the database.

  • sql_q (str) – Sql query to use. None will use a default query to get all star info.

  • dtypes (list [str, dtype]) – Data types expected from the database. None will use the defaults.

  • db_address (str) – Database data path. Default uses db in $RUBIN_SIM_DATA/skybrightness.

  • filt (str) – Filter in which to fetch stellar observation data.

Returns:

  • data (np.ndarray, (N,)) – Stellar observation data.

  • mjd (float) – MJD of the observations.

rubin_sim.skybrightness.diode_sky_db(mid_mjd, sql_q=None, dtypes=None, db_address=None, clean=True)[source]

Fetch diode measurements of skybrightness.

rubin_sim.skybrightness.id2intid(ids)[source]

Take an array of ids, and convert them to an integer id. Handy if you want to put things into a sparse array.

rubin_sim.skybrightness.intid2id(intids, uintids, uids, dtype=<class 'int'>)[source]

convert an int back to an id

rubin_sim.skybrightness.just_return(inval)[source]

Really, just return the input.

Parameters:

input (anything)

Returns:

input – Just return whatever you sent in.

Return type:

anything

rubin_sim.skybrightness.load_spec_files(filenames, mags=False)[source]

Load up the ESO spectra.

The ESO npz files contain the following arrays: * filter_wave: The central wavelengths of the pre-computed magnitudes * wave: wavelengths for the spectra * spec: array of spectra and magnitudes along with the relevant variable inputs. For example, airglow has dtype = [(‘airmass’, ‘<f8’), (‘solar_flux’, ‘<f8’), (‘spectra’,’<f8’, (17001,)), (‘mags’, ‘<f8’, (6,)]

For each unique airmass and solar_flux value, there is a 17001 elements spectra and 6 magnitudes.

rubin_sim.skybrightness.recalc_mags(data_dir=None)[source]

Recalculate the magnitudes for sky brightness components.

DANGER: Overwrites data files in place. The rubin_sim_data/skybrightness folder will need to be packaged and updated after running this to propagate changes to other users.

rubin_sim.skybrightness.robust_rms(array, missing=0.0)[source]

Use the interquartile range to compute a robust approximation of the RMS. if passed an array smaller than 2 elements, return missing value

rubin_sim.skybrightness.simple_twi(xdata, *args)[source]

Fit a simple slope and constant to many healpixels

xdata should have keys: sunAlt hpid

args: 0: slope 1:hpid: magnitudes hpid+1:2*hpid: constant offsets

rubin_sim.skybrightness.twilight_func(xdata, *args, amCut=1.0)[source]

xdata: numpy array with columns ‘alt’, ‘az’, ‘sunAlt’ all in radians. az should be relative to the sun (i.e., sun is at az zero.

based on what I’ve seen, here’s my guess for how to fit the twilight: args[0] = ratio of (zenith twilight flux at sun_alt = -12) and dark sky zenith flux args[1] = decay slope for all pixels (mags/radian) args[2] = airmass term for hemisphere away from the sun. (factor to multiply max brightness at zenith by) args[3] = az term for hemisphere towards sun args[4] = zenith dark sky flux args[5:] = zenith dark sky times constant (optional)

amCutfloat (1.0)

The airmass cut to apply to use only the away from sun fit. Was set to 1.1 previously for not very clear reasons.

rubin_sim.skybrightness.wrap_ra(ra)[source]

Wrap only RA values into 0-2pi (using mod).

rubin_sim.skybrightness.zenith_twilight(alpha, *args)[source]

The flux at zenith as a linear combination of a twilight component and a constant: alpha = sun altitude (radians) args[0] = ratio of (zenith twilight flux at sunAlt = -12) and dark sky zenith flux args[1] = decay slope for all pixels (mags/radian) args[2] = airmass term for hemisphere away from the sun. (factor to multiply max brightness at zenith by) args[3] = az term for hemisphere towards sun args[4] = zenith dark sky flux